Matrices

Matrices — construct and manipulate matrices

Functions

#define gretl_matrix_get()
#define gretl_cmatrix_get()
#define gretl_vector_get()
#define gretl_matrix_set()
#define gretl_cmatrix_set()
#define gretl_vector_set()
#define gretl_matrix_cols()
#define gretl_matrix_rows()
#define gretl_vector_get_length()
#define gretl_vector_alloc()
#define gretl_column_vector_alloc()
#define gretl_vector_free()
#define gretl_matrix_is_scalar()
#define gretl_matrix_is_cscalar()
#define gretl_is_null_matrix()
#define gretl_is_complex()
int gretl_matrix_set_complex ()
int gretl_matrix_set_complex_full ()
int get_gretl_matrix_err ()
void clear_gretl_matrix_err ()
void gretl_matrix_print ()
void gretl_matrix_print2 ()
int gretl_matrix_na_check ()
int gretl_matrix_is_symmetric ()
int gretl_matrix_is_idempotent ()
void gretl_matrix_xtr_symmetric ()
void gretl_matrix_set_equals_tolerance ()
void gretl_matrix_unset_equals_tolerance ()
gretl_matrix * gretl_matrix_alloc ()
gretl_matrix * gretl_cmatrix_new ()
gretl_matrix * gretl_cmatrix_new0 ()
gretl_matrix * gretl_cmatrix_new1 ()
gretl_matrix * gretl_matching_matrix_new ()
gretl_matrix * gretl_matrix_reuse ()
int gretl_matrix_realloc ()
gretl_matrix * gretl_matrix_init ()
gretl_matrix * gretl_matrix_init_full ()
gretl_matrix * gretl_matrix_replace ()
int gretl_matrix_replace_content ()
void gretl_matrix_block_destroy ()
void gretl_matrix_block_zero ()
gretl_matrix_block * gretl_matrix_block_new ()
int gretl_matrix_block_n_matrices ()
gretl_matrix * gretl_matrix_block_get_matrix ()
gretl_matrix * gretl_identity_matrix_new ()
gretl_matrix * gretl_DW_matrix_new ()
gretl_matrix * gretl_zero_matrix_new ()
gretl_matrix * gretl_unit_matrix_new ()
gretl_matrix * gretl_null_matrix_new ()
gretl_matrix * gretl_matrix_seq ()
gretl_matrix * gretl_matrix_copy ()
int gretl_matrix_copy_row ()
int gretl_matrix_inscribe_I ()
gretl_matrix * gretl_matrix_copy_transpose ()
gretl_matrix * gretl_matrix_reverse_rows ()
gretl_matrix * gretl_matrix_reverse_cols ()
gretl_matrix * gretl_matrix_get_diagonal ()
int gretl_matrix_set_diagonal ()
gretl_matrix * gretl_matrix_get_triangle ()
int gretl_matrix_set_triangle ()
int gretl_matrix_get_row ()
double gretl_matrix_trace ()
int gretl_matrix_random_fill ()
gretl_matrix * gretl_random_matrix_new ()
gretl_matrix * gretl_matrix_resample ()
int gretl_matrix_resample2 ()
gretl_matrix * gretl_matrix_block_resample ()
int gretl_matrix_block_resample2 ()
double gretl_vector_mean ()
double gretl_vector_variance ()
void gretl_matrix_zero ()
int gretl_matrix_zero_upper ()
int gretl_matrix_zero_lower ()
int gretl_matrix_mirror ()
void gretl_matrix_fill ()
void gretl_matrix_multiply_by_scalar ()
int gretl_matrix_divide_by_scalar ()
void gretl_matrix_switch_sign ()
gretl_matrix * gretl_matrix_dot_op ()
ConfType dot_operator_conf ()
gretl_matrix * gretl_matrix_complex_multiply ()
gretl_matrix * gretl_matrix_divide ()
gretl_matrix * gretl_matrix_complex_divide ()
gretl_matrix * gretl_matrix_exp ()
gretl_matrix * gretl_matrix_polroots ()
void gretl_matrix_raise ()
void gretl_matrix_free ()
double * gretl_matrix_steal_data ()
int gretl_vector_copy_values ()
int gretl_matrix_copy_values ()
int gretl_matrix_copy_data ()
int gretl_matrix_copy_values_shaped ()
int gretl_matrix_add_to ()
int gretl_matrix_add ()
int gretl_matrix_add_transpose_to ()
int gretl_matrix_subtract_from ()
int gretl_matrix_subtract ()
int gretl_matrix_subtract_reversed ()
int gretl_matrix_I_minus ()
int gretl_matrix_transpose_in_place ()
int gretl_matrix_transpose ()
int gretl_square_matrix_transpose ()
int gretl_matrix_add_self_transpose ()
int gretl_matrix_vectorize ()
gretl_matrix * gretl_matrix_vectorize_new ()
int gretl_matrix_unvectorize ()
int gretl_matrix_vectorize_h ()
int gretl_matrix_vectorize_h_skip ()
int gretl_matrix_unvectorize_h ()
int gretl_matrix_unvectorize_h_diag ()
int gretl_matrix_inscribe_matrix ()
int gretl_matrix_extract_matrix ()
int gretl_matrix_multiply_mod ()
int gretl_matrix_multiply_mod_single ()
int gretl_matrix_multiply ()
int gretl_matrix_multiply_single ()
gretl_matrix * gretl_matrix_multiply_new ()
void gretl_blas_dgemm ()
void gretl_blas_dsymm ()
int gretl_matrix_kronecker_product ()
gretl_matrix * gretl_matrix_kronecker_product_new ()
int gretl_matrix_hdproduct ()
gretl_matrix * gretl_matrix_hdproduct_new ()
int gretl_matrix_I_kronecker ()
gretl_matrix * gretl_matrix_I_kronecker_new ()
int gretl_matrix_kronecker_I ()
gretl_matrix * gretl_matrix_kronecker_I_new ()
gretl_matrix * gretl_matrix_pow ()
double gretl_matrix_dot_product ()
double gretl_vector_dot_product ()
gretl_matrix * gretl_rmatrix_vector_stat ()
gretl_matrix * gretl_matrix_column_sd ()
void gretl_matrix_demean_by_row ()
int gretl_matrix_standardize ()
int gretl_matrix_center ()
gretl_matrix * gretl_matrix_quantiles ()
double gretl_matrix_determinant ()
double gretl_matrix_log_determinant ()
double gretl_matrix_log_abs_determinant ()
double gretl_vcv_log_determinant ()
double gretl_matrix_one_norm ()
double gretl_matrix_infinity_norm ()
int gretl_LU_solve ()
int gretl_LU_solve_invert ()
int gretl_cholesky_decomp_solve ()
int gretl_cholesky_solve ()
int gretl_cholesky_invert ()
gretl_vector * gretl_toeplitz_solve ()
gretl_matrix * gretl_matrix_XTX_new ()
int gretl_inverse_from_cholesky_decomp ()
int gretl_invert_general_matrix ()
int gretl_invert_symmetric_indef_matrix ()
int gretl_invert_symmetric_matrix ()
int gretl_invert_symmetric_matrix2 ()
int gretl_invert_packed_symmetric_matrix ()
int gretl_invert_triangular_matrix ()
int gretl_invert_diagonal_matrix ()
int gretl_invert_matrix ()
int gretl_matrix_moore_penrose ()
int gretl_SVD_invert_matrix ()
int gretl_invpd ()
int gretl_matrix_SVD ()
double gretl_symmetric_matrix_rcond ()
double gretl_matrix_rcond ()
double gretl_matrix_cond_index ()
int gretl_symmetric_eigen_sort ()
gretl_matrix * gretl_general_matrix_eigenvals ()
gretl_matrix * gretl_symmetric_matrix_eigenvals ()
double gretl_symmetric_matrix_min_eigenvalue ()
gretl_matrix * gretl_symm_matrix_eigenvals_descending ()
gretl_matrix * gretl_gensymm_eigenvals ()
gretl_matrix * gretl_dgeev ()
gretl_matrix * old_eigengen ()
double gretl_symm_matrix_lambda_min ()
double gretl_symm_matrix_lambda_max ()
gretl_matrix * gretl_matrix_right_nullspace ()
gretl_matrix * gretl_matrix_left_nullspace ()
gretl_matrix * gretl_matrix_row_concat ()
gretl_matrix * gretl_matrix_col_concat ()
gretl_matrix * gretl_matrix_direct_sum ()
int gretl_matrix_inplace_colcat ()
gretl_matrix * gretl_matrix_cumcol ()
gretl_matrix * gretl_matrix_diffcol ()
gretl_matrix * gretl_matrix_lag ()
int gretl_matrix_inplace_lag ()
int gretl_matrix_cholesky_decomp ()
int gretl_matrix_psd_root ()
int gretl_matrix_QR_decomp ()
int gretl_matrix_QR_pivot_decomp ()
double gretl_triangular_matrix_rcond ()
int gretl_check_QR_rank ()
int gretl_matrix_rank ()
int gretl_matrix_ols ()
int gretl_matrix_multi_ols ()
int gretl_matrix_multi_SVD_ols ()
int gretl_matrix_QR_ols ()
double gretl_matrix_r_squared ()
int gretl_matrix_SVD_johansen_solve ()
int gretl_matrix_restricted_ols ()
int gretl_matrix_restricted_multi_ols ()
int gretl_matrix_SVD_ols ()
int gretl_matrix_qform ()
int gretl_matrix_diag_qform ()
double gretl_scalar_qform ()
int gretl_matrix_columnwise_product ()
int gretl_matrix_diagonal_sandwich ()
int gretl_matrix_set_t1 ()
int gretl_matrix_set_t2 ()
int gretl_matrix_get_t1 ()
int gretl_matrix_get_t2 ()
int gretl_matrix_is_dated ()
int gretl_is_identity_matrix ()
int gretl_is_zero_matrix ()
gretl_matrix * gretl_matrix_isfinite ()
int gretl_matrix_get_structure ()
int gretl_matrices_are_equal ()
gretl_matrix * gretl_covariance_matrix ()
gretl_matrix * gretl_matrix_GG_inverse ()
gretl_matrix * gretl_matrix_varsimul ()
gretl_matrix ** gretl_matrix_array_new ()
gretl_matrix ** gretl_matrix_array_new_with_size ()
void gretl_matrix_array_free ()
gretl_matrix * gretl_matrix_values ()
gretl_matrix * gretl_matrix_shape ()
gretl_matrix * gretl_matrix_trim_rows ()
gretl_matrix * gretl_matrix_minmax ()
double gretl_matrix_global_minmax ()
double gretl_matrix_global_sum ()
gretl_matrix * gretl_matrix_pca ()
gretl_matrix * gretl_matrix_xtab ()
gretl_matrix * matrix_matrix_xtab ()
gretl_matrix * gretl_matrix_bool_sel ()
gretl_matrix * gretl_matrix_sort_by_column ()
gretl_matrix * gretl_vector_sort ()
gretl_matrix * gretl_matrix_covariogram ()
gretl_matrix * gretl_matrix_commute ()
void gretl_matrix_transcribe_obs_info ()
int gretl_matrix_set_colnames ()
int gretl_matrix_set_rownames ()
const char ** gretl_matrix_get_colnames ()
const char ** gretl_matrix_get_rownames ()
void gretl_matrix_destroy_info ()
void lapack_mem_free ()
void set_blas_mnk_min ()
int get_blas_mnk_min ()
void set_simd_k_max ()
int get_simd_k_max ()
void set_simd_mn_min ()
int get_simd_mn_min ()

Types and Values

Includes

#include <libgretl.h>

Description

Libgretl implements most of the matrix functionality that is likely to be required in econometric calculation. For basics such as decomposition and inversion we use LAPACK as the underlying engine.

To get yourself a gretl matrix, use gretl_matrix_alloc() or one of the more specialized constructors; to free such a matrix use gretl_matrix_free().

Functions

gretl_matrix_get()

#define gretl_matrix_get(m,i,j) (m->val[(j)*m->rows+(i)])

Parameters

m

matrix.

 

i

row.

 

j

column.

 

Returns

the i , j element of m .


gretl_cmatrix_get()

#define gretl_cmatrix_get(m,i,j) (m->z[(j)*m->rows+(i)])

gretl_vector_get()

#define gretl_vector_get(v,i) (v->val[i])

Parameters

v

vector.

 

i

index.

 

Returns

element i of v .


gretl_matrix_set()

#define gretl_matrix_set(m,i,j,x) ((m)->val[(j)*(m)->rows+(i)]=x)

Sets the i , j element of m to x .

Parameters

m

matrix.

 

i

row.

 

j

column.

 

x

value to set.

 

gretl_cmatrix_set()

#define gretl_cmatrix_set(m,i,j,x) ((m)->z[(j)*(m)->rows+(i)]=x)

gretl_vector_set()

#define gretl_vector_set(v,i,x) ((v)->val[i]=x)

Sets element i of v to x .

Parameters

v

vector.

 

i

index.

 

x

value to set.

 

gretl_matrix_cols()

#define gretl_matrix_cols(m) ((m == NULL)? 0 : m->cols)

Parameters

m

matrix to query.

 

Returns

the number of columns in m .


gretl_matrix_rows()

#define gretl_matrix_rows(m) ((m == NULL)? 0 : m->rows)

Parameters

m

matrix to query.

 

Returns

the number of rows in m .


gretl_vector_get_length()

#define             gretl_vector_get_length(v)

Parameters

v

vector to examine.

 

Returns

the length of vector v (without regard to whether it is a row or column vector).


gretl_vector_alloc()

#define gretl_vector_alloc(i) gretl_matrix_alloc(1,(i))

Parameters

i

number of columns.

 

Returns

a new gretl_vector with i columns.


gretl_column_vector_alloc()

#define gretl_column_vector_alloc(i) gretl_matrix_alloc((i),1)

Parameters

i

number of rows.

 

Returns

a new column gretl_vector with i rows.


gretl_vector_free()

#define gretl_vector_free(v) gretl_matrix_free(v)

Frees the vector v and its associated storage.

Parameters

v

gretl_vector to free.

 

gretl_matrix_is_scalar()

#define             gretl_matrix_is_scalar(m)

Parameters

m

matrix to test.

 

Returns

1 if m is 1 x 1, else 0.


gretl_matrix_is_cscalar()

#define             gretl_matrix_is_cscalar(m)

gretl_is_null_matrix()

#define gretl_is_null_matrix(m) (m == NULL || m->rows == 0 || m->cols == 0)

gretl_is_complex()

#define gretl_is_complex(m) (m != NULL && m->is_complex == 1)

gretl_matrix_set_complex ()

int
gretl_matrix_set_complex (gretl_matrix *m,
                          int c);

gretl_matrix_set_complex_full ()

int
gretl_matrix_set_complex_full (gretl_matrix *m,
                               int c);

get_gretl_matrix_err ()

int
get_gretl_matrix_err (void);

clear_gretl_matrix_err ()

void
clear_gretl_matrix_err (void);

gretl_matrix_print ()

void
gretl_matrix_print (const gretl_matrix *m,
                    const char *msg);

Prints the given matrix to stderr.

Parameters

m

matrix.

 

msg

message to print with matrix, or NULL.

 

gretl_matrix_print2 ()

void
gretl_matrix_print2 (const gretl_matrix *m,
                     const char *msg);

Prints the given matrix to stdout.

Parameters

m

matrix.

 

msg

message to print with matrix, or NULL.

 

gretl_matrix_na_check ()

int
gretl_matrix_na_check (const gretl_matrix *m);

gretl_matrix_is_symmetric ()

int
gretl_matrix_is_symmetric (const gretl_matrix *m);

Parameters

m

gretl_matrix.

 

Returns

1 if m is symmetric (with a small relative tolerance for asymmetry), otherwise 0.


gretl_matrix_is_idempotent ()

int
gretl_matrix_is_idempotent (const gretl_matrix *m,
                            double tol);

Parameters

m

gretl_matrix.

 

tol

numerical tolerance

 

Returns

1 if m is idempotent, otherwise 0.


gretl_matrix_xtr_symmetric ()

void
gretl_matrix_xtr_symmetric (gretl_matrix *m);

Computes the symmetric part of m by averaging its off-diagonal elements.

Parameters

m

gretl_matrix.

 

gretl_matrix_set_equals_tolerance ()

void
gretl_matrix_set_equals_tolerance (double tol);

Sets the tolerance for judging whether or not a matrix is symmetric (see gretl_matrix_is_symmetric() and also gretl_invert_symmetric_matrix()). The tolerance is the maximum relative difference between corresponding off-diagonal elements that is acceptable in a supposedly "symmetric" matrix. The default value is 1.0e-9.

Parameters

tol

tolerance value.

 

gretl_matrix_unset_equals_tolerance ()

void
gretl_matrix_unset_equals_tolerance (void);

Sets the tolerance for judging whether or not a matrix is symmetric to its default value. See also gretl_matrix_set_equals_tolerance().


gretl_matrix_alloc ()

gretl_matrix *
gretl_matrix_alloc (int rows,
                    int cols);

Parameters

rows

desired number of rows in matrix.

 

cols

desired number of columns.

 

Returns

pointer to a newly allocated gretl_matrix, or NULL on failure. Note that the actual data storage is not initialized.


gretl_cmatrix_new ()

gretl_matrix *
gretl_cmatrix_new (int r,
                   int c);

gretl_cmatrix_new0 ()

gretl_matrix *
gretl_cmatrix_new0 (int r,
                    int c);

gretl_cmatrix_new1 ()

gretl_matrix *
gretl_cmatrix_new1 (int r,
                    int c);

gretl_matching_matrix_new ()

gretl_matrix *
gretl_matching_matrix_new (int r,
                           int c,
                           const gretl_matrix *m);

gretl_matrix_reuse ()

gretl_matrix *
gretl_matrix_reuse (gretl_matrix *m,
                    int rows,
                    int cols);

An "experts only" memory-conservation trick. If m is an already-allocated gretl matrix, you can "resize" it by specifying a new number of rows and columns. This works only if the product of rows and cols is less than or equal to the product of the number of rows and columns in the matrix as originally allocated; no actual reallocation of memory is performed. If you "reuse" with an excessive number of rows or columns you will surely crash your program or smash the stack. Note also that the matrix-pointer returned is not really new, and when the matrix is to be freed, gretl_matrix_free() should be applied only once.

Parameters

m

matrix to reuse.

 

rows

desired number of rows in "new" matrix, or -1 to leave the current value unchanged.

 

cols

desired number of columns in "new" matrix, or -1 to leave the current value unchanged.

 

Returns

pointer to the "resized" gretl_matrix.


gretl_matrix_realloc ()

int
gretl_matrix_realloc (gretl_matrix *m,
                      int rows,
                      int cols);

Reallocates the storage in m to the specified dimensions.

Parameters

m

matrix to reallocate.

 

rows

desired number of rows in "new" matrix.

 

cols

desired number of columns in "new" matrix.

 

Returns

0 on success, E_ALLOC on failure.


gretl_matrix_init ()

gretl_matrix *
gretl_matrix_init (gretl_matrix *m);

Initializes m to be zero by zero with NULL data.

Parameters

m

matrix to be initialized.

 

gretl_matrix_init_full ()

gretl_matrix *
gretl_matrix_init_full (gretl_matrix *m,
                        int rows,
                        int cols,
                        double *val);

Initializes m to be rows by cols and have data val . This intended for use with automatic matrices declared "on the stack". It is up to the user to ensure that the size of val is compatible with the rows and cols specification.

Parameters

m

matrix to be initialized.

 

rows

number of rows.

 

cols

number of columns.

 

val

data array.

 

gretl_matrix_replace ()

gretl_matrix *
gretl_matrix_replace (gretl_matrix **pa,
                      gretl_matrix *b);

Frees the matrix at location pa and substitutes b .

Parameters

pa

location of matrix to be replaced.

 

b

replacement matrix.

 

Returns

the replacement matrix.


gretl_matrix_replace_content ()

int
gretl_matrix_replace_content (gretl_matrix *targ,
                              gretl_matrix *donor);

Moves the content of donor into targ ; donor becomes a null matrix in consequence.

Parameters

targ

matrix to receive new content.

 

donor

matrix to donate content.

 

Returns

0 on success, non-zero on error.


gretl_matrix_block_destroy ()

void
gretl_matrix_block_destroy (gretl_matrix_block *B);

gretl_matrix_block_zero ()

void
gretl_matrix_block_zero (gretl_matrix_block *B);

gretl_matrix_block_new ()

gretl_matrix_block *
gretl_matrix_block_new (gretl_matrix **pm,
                        ...);

gretl_matrix_block_n_matrices ()

int
gretl_matrix_block_n_matrices (gretl_matrix_block *B);

gretl_matrix_block_get_matrix ()

gretl_matrix *
gretl_matrix_block_get_matrix (gretl_matrix_block *B,
                               int i);

gretl_identity_matrix_new ()

gretl_matrix *
gretl_identity_matrix_new (int n);

Parameters

n

desired number of rows and columns in the matrix.

 

Returns

pointer to a newly allocated identity matrix, or NULL on failure.


gretl_DW_matrix_new ()

gretl_matrix *
gretl_DW_matrix_new (int n);

Parameters

n

desired number of rows and columns in the matrix.

 

Returns

pointer to a newly allocated Durbin-Watson matrix, or NULL on failure. This is a tridiagonal matrix with 2 on the leading diagonal (apart from 1 at the ends) and -1 on the supra- and infra-diagonals.


gretl_zero_matrix_new ()

gretl_matrix *
gretl_zero_matrix_new (int r,
                       int c);

Parameters

r

desired number of rows in the matrix.

 

c

desired number of columns in the matrix.

 

Returns

pointer to a newly allocated zero matrix, or NULL on failure.


gretl_unit_matrix_new ()

gretl_matrix *
gretl_unit_matrix_new (int r,
                       int c);

Parameters

r

desired number of rows in the matrix.

 

c

desired number of columns in the matrix.

 

Returns

pointer to a newly allocated matrix, all of whose elements equal 1, or NULL on failure.


gretl_null_matrix_new ()

gretl_matrix *
gretl_null_matrix_new (void);

Returns

pointer to a newly allocated null matrix, or NULL on failure.


gretl_matrix_seq ()

gretl_matrix *
gretl_matrix_seq (double start,
                  double end,
                  double step,
                  int *err);

Parameters

start

first element.

 

end

last element.

 

step

positive step.

 

err

location to recieve error code.

 

Returns

pointer to a row vector, containing values from start to end , in decreasing order if start > end -- or NULL on failure.


gretl_matrix_copy ()

gretl_matrix *
gretl_matrix_copy (const gretl_matrix *m);

Parameters

m

source matrix to be copied.

 

Returns

an allocated copy of matrix m , or NULL on failure.


gretl_matrix_copy_row ()

int
gretl_matrix_copy_row (gretl_matrix *dest,
                       int di,
                       const gretl_matrix *src,
                       int si);

Copies the values from row si of src into row di of dest , provided src and dest have the same number of columns.

Parameters

dest

destination matrix.

 

di

row to copy into.

 

src

source matrix.

 

si

row to copy from.

 

Returns

0 on success, non-zero on failure.


gretl_matrix_inscribe_I ()

int
gretl_matrix_inscribe_I (gretl_matrix *m,
                         int row,
                         int col,
                         int n);

Writes an n x n identity matrix into matrix m , the top left-hand corner of the insertion being given by row and col (which are 0-based).

Parameters

m

original matrix.

 

row

top row for insertion.

 

col

leftmost row for insertion.

 

n

dimension (rows and columns) of identity matrix.

 

Returns

0 on successful completion, or E_NONCONF if an identity matrix of the specified size cannot be fitted into m at the specified location.


gretl_matrix_copy_transpose ()

gretl_matrix *
gretl_matrix_copy_transpose (const gretl_matrix *m);

Parameters

m

source matrix to be copied.

 

Returns

an allocated copy of the tranpose of m , or NULL on failure.


gretl_matrix_reverse_rows ()

gretl_matrix *
gretl_matrix_reverse_rows (const gretl_matrix *m,
                           int *err);

Parameters

m

source matrix whose rows are to be reversed.

 

err

location to receive error code.

 

Returns

a matrix with the same rows as m , last to first.


gretl_matrix_reverse_cols ()

gretl_matrix *
gretl_matrix_reverse_cols (const gretl_matrix *m,
                           int *err);

Parameters

m

source matrix whose columns are to be reversed.

 

err

location to receive error code.

 

Returns

a matrix with the same columns as m , last to first.


gretl_matrix_get_diagonal ()

gretl_matrix *
gretl_matrix_get_diagonal (const gretl_matrix *m,
                           int *err);

Parameters

m

input matrix.

 

err

location to receive error code.

 

Returns

a column vector containing the diagonal elements of m , otherwise NULL. A non-zero value is assigned via err on failure.


gretl_matrix_set_diagonal ()

int
gretl_matrix_set_diagonal (gretl_matrix *targ,
                           const gretl_matrix *src,
                           double x);

Sets the diagonal elements of targ using the elements of src , if non-NULL, or otherwise the constant value x . If src is given it must be a vector of length equal to that of the diagonal of targ (that is, the minimum of its rows and columns).

Parameters

targ

target matrix.

 

src

source vector (or NULL).

 

x

(alternative) source scalar.

 

Returns

0 on success, error code on non-conformability.


gretl_matrix_get_triangle ()

gretl_matrix *
gretl_matrix_get_triangle (const gretl_matrix *m,
                           int upper,
                           int *err);

Parameters

m

source matrix (real or complex).

 

upper

flag to get the upper part, the default being to get the lower part.

 

err

location to receive error code.

 

Returns

A column vector holding the vec of either the infra- or supra-diagonal elements of m , or NULL on failure. Note that the "part" returned may not be an actual triangle if m is not square.


gretl_matrix_set_triangle ()

int
gretl_matrix_set_triangle (gretl_matrix *targ,
                           const gretl_matrix *src,
                           double x,
                           int upper);

Sets the lower or upper elements of the matrix targ using the elements of src , if non-NULL, or otherwise the constant value x .

If src is given it must be a vector of length equal to that of the number of infra- or supra-diagonal elements of targ .

Parameters

targ

target matrix.

 

src

source vector (or NULL).

 

x

(alternative) source scalar.

 

upper

flag to set the upper part, the default being to set the lower.

 

Returns

0 on success, error code on non-conformability.


gretl_matrix_get_row ()

int
gretl_matrix_get_row (const gretl_matrix *m,
                      int i,
                      gretl_vector *v);

Copies row i of matrix m into vector v .

Parameters

m

input matrix.

 

i

index of row to access.

 

v

location to receive row values.

 

Returns

0 on success, non-zero on error.


gretl_matrix_trace ()

double
gretl_matrix_trace (const gretl_matrix *m);

Parameters

m

square input matrix.

 

Returns

the trace (sum of diagonal elements) of m , if m is square, otherwise NADBL.


gretl_matrix_random_fill ()

int
gretl_matrix_random_fill (gretl_matrix *m,
                          int dist);

Fills m with pseudo-random values from either the uniform or the standard normal distribution.

Parameters

m

input matrix.

 

dist

either D_UNIFORM or D_NORMAL.

 

Returns

0 on success, 1 on failure.


gretl_random_matrix_new ()

gretl_matrix *
gretl_random_matrix_new (int r,
                         int c,
                         int dist);

Creates a new $r x c matrix and filles it with pseudo-random values from either the uniform or the standard normal distribution.

Parameters

r

number of rows.

 

c

number of columns.

 

dist

either D_UNIFORM or D_NORMAL.

 

Returns

allocated matrix or NULL on failure.


gretl_matrix_resample ()

gretl_matrix *
gretl_matrix_resample (const gretl_matrix *m,
                       int draws,
                       int *err);

Parameters

m

input matrix.

 

draws

number of draws (or 0 to use the number or rows in m ).

 

err

location to receive error code.

 

Returns

a new matrix consisting of a random re-sampling (with replacement) of the rows of m , or NULL on failure.


gretl_matrix_resample2 ()

int
gretl_matrix_resample2 (gretl_matrix *targ,
                        const gretl_matrix *src);

gretl_matrix_block_resample ()

gretl_matrix *
gretl_matrix_block_resample (const gretl_matrix *m,
                             int blocklen,
                             int draws,
                             int *err);

Parameters

m

input matrix.

 

blocklen

length of moving blocks.

 

draws

number of draws (or 0 to use the rows of m ).

 

err

location to receive error code.

 

Returns

a new matrix consisting of a random re-sampling (with replacement) of the rows of m , using blocks of contiguous rows of length blocklen , or NULL on failure.


gretl_matrix_block_resample2 ()

int
gretl_matrix_block_resample2 (gretl_matrix *targ,
                              const gretl_matrix *src,
                              int blocklen,
                              int *z);

An "in-place" version of gretl_matrix_block_resample(). It is assumed that targ is a matrix of the same dimensions as src , that blocklen is greater than 1, and that z is long enough to hold n integers, where n is the number of rows in src divided by blocklen , rounded up to the nearest integer.

Parameters

targ

target matrix.

 

src

source matrix.

 

blocklen

length of moving blocks.

 

z

array of length XXX.

 

Returns

0 on success, non-zero on failure.


gretl_vector_mean ()

double
gretl_vector_mean (const gretl_vector *v);

Parameters

v

input vector.

 

Returns

the arithmetic mean of the elements of v , or NADBL on failure.


gretl_vector_variance ()

double
gretl_vector_variance (const gretl_vector *v);

Parameters

v

input vector.

 

Returns

the variance of the elements of v , or NADBL on failure.


gretl_matrix_zero ()

void
gretl_matrix_zero (gretl_matrix *m);

Sets all elements of m to zero.

Parameters

m

matrix to be set to zero.

 

gretl_matrix_zero_upper ()

int
gretl_matrix_zero_upper (gretl_matrix *m);

Sets the elements of m outside of the lower triangle to zero.

Parameters

m

square matrix to operate on.

 

Returns

0 on success, non-zero error code otherwise.


gretl_matrix_zero_lower ()

int
gretl_matrix_zero_lower (gretl_matrix *m);

Sets the elements of m outside of the upper triangle to zero.

Parameters

m

square matrix to operate on.

 

Returns

0 on success, non-zero error code otherwise.


gretl_matrix_mirror ()

int
gretl_matrix_mirror (gretl_matrix *m,
                     char uplo);

If uplo = 'L', copy the lower triangle of m into the upper triangle; or if uplo = 'U' copy the upper triangle into the lower, in either case producing a symmetric result.

Parameters

m

matrix to expand.

 

uplo

'L' or 'U'.

 

Returns

0 on success; non-zero error code if m is not square.


gretl_matrix_fill ()

void
gretl_matrix_fill (gretl_matrix *m,
                   double x);

Sets all entries in m to the value x .

Parameters

m

matrix to fill.

 

x

value with which to fill.

 

gretl_matrix_multiply_by_scalar ()

void
gretl_matrix_multiply_by_scalar (gretl_matrix *m,
                                 double x);

Multiplies all elements of m by x .

Parameters

m

matrix to operate on.

 

x

scalar by which to multiply.

 

gretl_matrix_divide_by_scalar ()

int
gretl_matrix_divide_by_scalar (gretl_matrix *m,
                               double x);

Divides all elements of m by x .

Parameters

m

matrix to operate on.

 

x

scalar by which to divide.

 

Returns

0 on success, 1 if x = 0.


gretl_matrix_switch_sign ()

void
gretl_matrix_switch_sign (gretl_matrix *m);

Changes the sign of each element of m .

Parameters

m

matrix to operate on.

 

gretl_matrix_dot_op ()

gretl_matrix *
gretl_matrix_dot_op (const gretl_matrix *a,
                     const gretl_matrix *b,
                     int op,
                     int *err);

Parameters

a

left-hand matrix.

 

b

right-hand matrix.

 

op

operator.

 

err

location to receive error code.

 

Returns

a new matrix, each of whose elements is the result of (x op y), where x and y are the corresponding elements of the matrices a and b (or NULL on failure).


dot_operator_conf ()

ConfType
dot_operator_conf (const gretl_matrix *A,
                   const gretl_matrix *B,
                   int *r,
                   int *c);

Used to establish the dimensions of the result of a "dot" operation such as A .* B or A .+ B.

Parameters

A

first matrix.

 

B

second matrix.

 

r

pointer to rows of result.

 

c

pointer to columns of result.

 

Returns

a numeric code identifying the convention to be used; CONF_NONE indicates non-conformability.


gretl_matrix_complex_multiply ()

gretl_matrix *
gretl_matrix_complex_multiply (const gretl_matrix *a,
                               const gretl_matrix *b,
                               int force_complex,
                               int *err);

Computes the complex product of a and b . The first column in these matrices is assumed to contain real values, and the second column (if present) imaginary coefficients.

Parameters

a

m x (1 or 2) matrix.

 

b

m x (1 or 2) matrix.

 

force_complex

see below.

 

err

location to receive error code.

 

Returns

a matrix with the result of the multiplication of the two vectors of complex numbers. If both a and b have no imaginary part and the force_complex flag is zero, the return value will be m x 1, otherwise it will be m x 2.


gretl_matrix_divide ()

gretl_matrix *
gretl_matrix_divide (const gretl_matrix *a,
                     const gretl_matrix *b,
                     GretlMatrixMod mod,
                     int *err);

Follows the semantics of Matlab/Octave for left and right matrix "division". In left division, A \ B is in principle equivalent to A^{-1} * B, and in right division A / B is in principle equivalent to A * B^{-1}, but the result is obtained without explicit computation of the inverse.

In left division a and b must have the same number of rows; in right division they must have the same number of columns.

Parameters

a

left-hand matrix.

 

b

right-hand matrix.

 

mod

GRETL_MOD_NONE for left division, or GRETL_MOD_TRANSPOSE for right division.

 

err

location to receive error code.

 

Returns

the "quotient" matrix, or NULL on failure.


gretl_matrix_complex_divide ()

gretl_matrix *
gretl_matrix_complex_divide (const gretl_matrix *a,
                             const gretl_matrix *b,
                             int force_complex,
                             int *err);

Computes the complex division of a over b . The first column in these matrices is assumed to contain real values, and the second column (if present) imaginary coefficients.

Parameters

a

m x (1 or 2) matrix.

 

b

m x (1 or 2) matrix.

 

force_complex

see below.

 

err

location to receive error code.

 

Returns

a matrix with the result of the division of the two vectors of complex numbers. If both a and b have no imaginary part and the force_complex flag is zero, the return value will be m x 1, otherwise it will be m x 2.


gretl_matrix_exp ()

gretl_matrix *
gretl_matrix_exp (const gretl_matrix *m,
                  int *err);

Calculates the matrix exponential of m , using algorithm 11.3.1 from Golub and Van Loan, "Matrix Computations", 3e.

Parameters

m

square matrix to operate on.

 

err

location to receive error code.

 

Returns

the exponential, or NULL on failure.


gretl_matrix_polroots ()

gretl_matrix *
gretl_matrix_polroots (const gretl_matrix *a,
                       int force_complex,
                       int legacy,
                       int *err);

Calculates the roots of the polynomial with coefficients given by a . If the degree of the polynomial is p, then a should contain p + 1 coefficients in ascending order, i.e. starting with the constant and ending with the coefficient on x^p.

As a transitional measure, if legacy is non-zero, in the case of complex roots the return vector is in old-style format: real values in column 0 and imaginary parts in column 1.

Parameters

a

vector of coefficients.

 

force_complex

see below.

 

legacy

see below.

 

err

location to receive error code.

 

Returns

by default, a regular p-vector if all the roots are real, otherwise a p x 1 complex matrix. The force_complex flag can be used to produce a complex return value even if the imaginary parts are all zero.


gretl_matrix_raise ()

void
gretl_matrix_raise (gretl_matrix *m,
                    double x);

Raises each element of m to the power x .

Parameters

m

matrix to operate on.

 

x

exponent.

 

gretl_matrix_free ()

void
gretl_matrix_free (gretl_matrix *m);

Frees the allocated storage in m , then m itself.

Parameters

m

matrix to be freed.

 

gretl_matrix_steal_data ()

double *
gretl_matrix_steal_data (gretl_matrix *m);

"Steals" the allocated data from m , which is left with a NULL data pointer.

Parameters

m

matrix to operate on.

 

Returns

a pointer to the "stolen" data.


gretl_vector_copy_values ()

int
gretl_vector_copy_values (gretl_vector *targ,
                          const gretl_vector *src);

Copies the elements of src into the corresponding elements of targ .

Parameters

targ

target vector.

 

src

source vector.

 

Returns

0 on successful completion, or E_NONCONF if the two vectors are not of the same length.


gretl_matrix_copy_values ()

int
gretl_matrix_copy_values (gretl_matrix *targ,
                          const gretl_matrix *src);

Copies the elements of src into the corresponding elements of targ .

Parameters

targ

target matrix.

 

src

source matrix.

 

Returns

0 on successful completion, or E_NONCONF if the two matrices are not conformable for the operation.


gretl_matrix_copy_data ()

int
gretl_matrix_copy_data (gretl_matrix *targ,
                        const gretl_matrix *src);

Copies all data from src to targ : this does the same as gretl_matrix_copy_values() but also transcribes t1, t2, colnames and rownames information, if present.

Parameters

targ

target matrix.

 

src

source matrix.

 

Returns

0 on successful completion, non-zero code on failure.


gretl_matrix_copy_values_shaped ()

int
gretl_matrix_copy_values_shaped (gretl_matrix *targ,
                                 const gretl_matrix *src);

Copies the elements of src into targ , column by column.

Parameters

targ

target matrix.

 

src

source matrix.

 

Returns

0 on successful completion, or E_NONCONF if the two matrices do not contain the same number of elements.


gretl_matrix_add_to ()

int
gretl_matrix_add_to (gretl_matrix *targ,
                     const gretl_matrix *src);

Adds the elements of src to the corresponding elements of targ .

Parameters

targ

target matrix.

 

src

source matrix.

 

Returns

0 on successful completion, or E_NONCONF if the two matrices are not conformable for the operation. In the special case where src is in fact a scalar, the operation always goes through OK, with the scalar being added to each element of targ .


gretl_matrix_add ()

int
gretl_matrix_add (const gretl_matrix *a,
                  const gretl_matrix *b,
                  gretl_matrix *c);

Adds the elements of a and b , the result going to c .

Parameters

a

source matrix.

 

b

source matrix.

 

c

target matrix.

 

Returns

0 on successful completion, or E_NONCONF if the matrices are not conformable for the operation.


gretl_matrix_add_transpose_to ()

int
gretl_matrix_add_transpose_to (gretl_matrix *targ,
                               const gretl_matrix *src);

Adds the elements of src , transposed, to the corresponding elements of targ .

Parameters

targ

target matrix.

 

src

source matrix.

 

Returns

0 on successful completion, or E_NONCONF if the two matrices are not conformable for the operation.


gretl_matrix_subtract_from ()

int
gretl_matrix_subtract_from (gretl_matrix *targ,
                            const gretl_matrix *src);

Subtracts the elements of src from the corresponding elements of targ .

Parameters

targ

target matrix.

 

src

source matrix.

 

Returns

0 on successful completion, or E_NONCONF if the two matrices are not conformable for the operation. In the special case where src is in fact a scalar, the operation always goes through OK, with the scalar being subtracted from each element of targ .


gretl_matrix_subtract ()

int
gretl_matrix_subtract (const gretl_matrix *a,
                       const gretl_matrix *b,
                       gretl_matrix *c);

Subtracts the elements of b from the corresponding elements of a , the result going to c .

Parameters

a

source_matrix.

 

b

source matrix.

 

c

target matrix.

 

Returns

0 on successful completion, or E_NONCONF if the matrices are not conformable for the operation.


gretl_matrix_subtract_reversed ()

int
gretl_matrix_subtract_reversed (const gretl_matrix *a,
                                gretl_matrix *b);

Operates on b such that b_{ij} = a_{ij} - b_{ij}.

Parameters

a

m x n matrix.

 

b

m x n matrix.

 

Returns

0 on successful completion, or E_NONCONF if the two matrices are not conformable for the operation.


gretl_matrix_I_minus ()

int
gretl_matrix_I_minus (gretl_matrix *m);

Rewrites m as (I - m), where I is the n x n identity matrix.

Parameters

m

original square matrix, n x n.

 

Returns

0 on successful completion, or E_NONCONF if m is not square.


gretl_matrix_transpose_in_place ()

int
gretl_matrix_transpose_in_place (gretl_matrix *m);

Tranposes m in place.

Parameters

m

matrix to transpose.

 

Returns

0 on success, non-zero error code otherwise.


gretl_matrix_transpose ()

int
gretl_matrix_transpose (gretl_matrix *targ,
                        const gretl_matrix *src);

Fills out targ (which must be pre-allocated and of the right dimensions) with the transpose of src .

Parameters

targ

target matrix.

 

src

source matrix.

 

Returns

0 on success, non-zero error code otherwise.


gretl_square_matrix_transpose ()

int
gretl_square_matrix_transpose (gretl_matrix *m);

Transposes the matrix m .

Parameters

m

square matrix to operate on.

 

Returns

0 on success, non-zero error code otherwise.


gretl_matrix_add_self_transpose ()

int
gretl_matrix_add_self_transpose (gretl_matrix *m);

Adds the transpose of m to m itself, yielding a symmetric matrix.

Parameters

m

(square) matrix to operate on.

 

Returns

0 on successful completion, or 1 if the source matrix is not square.


gretl_matrix_vectorize ()

int
gretl_matrix_vectorize (gretl_matrix *targ,
                        const gretl_matrix *src);

Writes into targ vec(src ), that is, a column vector formed by stacking the columns of src .

Parameters

targ

target vector, (m * n) x 1.

 

src

source matrix, m x n.

 

Returns

0 on successful completion, or E_NONCONF if targ is not correctly dimensioned.


gretl_matrix_vectorize_new ()

gretl_matrix *
gretl_matrix_vectorize_new (const gretl_matrix *m);

Parameters

m

matrix to be vectorized.

 

Returns

a gretl column vector, vec(m ), or NULL on failure.


gretl_matrix_unvectorize ()

int
gretl_matrix_unvectorize (gretl_matrix *targ,
                          const gretl_matrix *src);

Writes successive blocks of length m from src into the successive columns of targ (that is, performs the inverse of the vec() operation).

Parameters

targ

target matrix, m x n.

 

src

source vector, (m * n) x 1.

 

Returns

0 on successful completion, or E_NONCONF if targ is not correctly dimensioned.


gretl_matrix_vectorize_h ()

int
gretl_matrix_vectorize_h (gretl_matrix *targ,
                          const gretl_matrix *src);

Writes into targ vech(src ), that is, a vector containing the lower-triangular elements of src . This is only useful for symmetric matrices, but for the sake of performance we don't check for that.

Parameters

targ

target vector, length n * (n+1)/2.

 

src

source square matrix, n x n.

 

Returns

0 on successful completion, or E_NONCONF if targ is not correctly dimensioned.


gretl_matrix_vectorize_h_skip ()

int
gretl_matrix_vectorize_h_skip (gretl_matrix *targ,
                               const gretl_matrix *src);

Like gretl_matrix_vectorize_h() except that the diagonal of src is skipped in the vectorization.

Parameters

targ

target vector, length n * (n-1)/2.

 

src

source square matrix, n x n.

 

Returns

0 on successful completion, or E_NONCONF if targ is not correctly dimensioned.


gretl_matrix_unvectorize_h ()

int
gretl_matrix_unvectorize_h (gretl_matrix *targ,
                            const gretl_matrix *src);

Rearranges successive blocks of decreasing length from src into the successive columns of targ (that is, performs the inverse of the vech() operation): targ comes out symmetric.

Parameters

targ

target matrix, n x n.

 

src

source vector, m x 1.

 

Returns

0 on successful completion, or E_NONCONF if targ is not correctly dimensioned.


gretl_matrix_unvectorize_h_diag ()

int
gretl_matrix_unvectorize_h_diag (gretl_matrix *targ,
                                 const gretl_matrix *src,
                                 double diag);

Like gretl_matrix_unvectorize_h(), but expects input src in which the diagonal elements are omitted, while diag gives the value to put in place of the omitted elements.

Parameters

targ

target matrix, n x n.

 

src

source vector.

 

diag

value for filling the diagonal.

 

Returns

0 on successful completion, or E_NONCONF if targ is not correctly dimensioned.


gretl_matrix_inscribe_matrix ()

int
gretl_matrix_inscribe_matrix (gretl_matrix *targ,
                              const gretl_matrix *src,
                              int row,
                              int col,
                              GretlMatrixMod mod);

Writes src into targ , starting at offset row , col . The targ matrix must be large enough to sustain the inscription of src at the specified point. If mod is GRETL_MOD_TRANSPOSE it is in fact the transpose of src that is written into targ .

Parameters

targ

target matrix.

 

src

source matrix.

 

row

row offset for insertion (0-based).

 

col

column offset for insertion.

 

mod

either GRETL_MOD_TRANSPOSE or GRETL_MOD_NONE.

 

Returns

0 on success, E_NONCONF if the matrices are not conformable for the operation.


gretl_matrix_extract_matrix ()

int
gretl_matrix_extract_matrix (gretl_matrix *targ,
                             const gretl_matrix *src,
                             int row,
                             int col,
                             GretlMatrixMod mod);

Writes into targ a sub-matrix of src , taken from the offset row , col . The targ matrix must be large enough to provide a sub-matrix of the dimensions of src . If mod is GRETL_MOD_TRANSPOSE it is in fact the transpose of the sub-matrix that that is written into targ .

Parameters

targ

target matrix.

 

src

source matrix.

 

row

row offset for extraction (0-based).

 

col

column offset for extraction.

 

mod

either GRETL_MOD_TRANSPOSE or GRETL_MOD_NONE.

 

Returns

0 on success, E_NONCONF if the matrices are not conformable for the operation.


gretl_matrix_multiply_mod ()

int
gretl_matrix_multiply_mod (const gretl_matrix *a,
                           GretlMatrixMod amod,
                           const gretl_matrix *b,
                           GretlMatrixMod bmod,
                           gretl_matrix *c,
                           GretlMatrixMod cmod);

Multiplies a (or a-transpose) into b (or b transpose), with the result written into c .

Parameters

a

left-hand matrix.

 

amod

modifier: GRETL_MOD_NONE or GRETL_MOD_TRANSPOSE.

 

b

right-hand matrix.

 

bmod

modifier: GRETL_MOD_NONE or GRETL_MOD_TRANSPOSE.

 

c

matrix to hold the product.

 

cmod

modifier: GRETL_MOD_NONE or GRETL_MOD_CUMULATE to add the result to the existing value of c , or GRETL_MOD_DECREMENT to subtract from the existing value of c .

 

Returns

0 on success; non-zero error code on failure.


gretl_matrix_multiply_mod_single ()

int
gretl_matrix_multiply_mod_single (const gretl_matrix *a,
                                  GretlMatrixMod amod,
                                  const gretl_matrix *b,
                                  GretlMatrixMod bmod,
                                  gretl_matrix *c,
                                  GretlMatrixMod cmod);

gretl_matrix_multiply ()

int
gretl_matrix_multiply (const gretl_matrix *a,
                       const gretl_matrix *b,
                       gretl_matrix *c);

Multiplies a into b , with the result written into c .

Parameters

a

left-hand matrix.

 

b

right-hand matrix.

 

c

matrix to hold the product.

 

Returns

0 on success; non-zero error code on failure.


gretl_matrix_multiply_single ()

int
gretl_matrix_multiply_single (const gretl_matrix *a,
                              const gretl_matrix *b,
                              gretl_matrix *c);

gretl_matrix_multiply_new ()

gretl_matrix *
gretl_matrix_multiply_new (const gretl_matrix *a,
                           const gretl_matrix *b,
                           int *err);

Multiplies a into b , with the result written into a newly allocated matrix.

Parameters

a

left-hand matrix.

 

b

right-hand matrix.

 

err

location for error code.

 

Returns

matrix product on success, or NULL on failure.


gretl_blas_dgemm ()

void
gretl_blas_dgemm (const gretl_matrix *a,
                  int atr,
                  const gretl_matrix *b,
                  int btr,
                  gretl_matrix *c,
                  GretlMatrixMod cmod,
                  int m,
                  int n,
                  int k);

gretl_blas_dsymm ()

void
gretl_blas_dsymm (const gretl_matrix *a,
                  int asecond,
                  const gretl_matrix *b,
                  int upper,
                  gretl_matrix *c,
                  GretlMatrixMod cmod,
                  int m,
                  int n);

gretl_matrix_kronecker_product ()

int
gretl_matrix_kronecker_product (const gretl_matrix *A,
                                const gretl_matrix *B,
                                gretl_matrix *K);

Writes the Kronecker product of A and B into K .

Parameters

A

left-hand matrix, p x q.

 

B

right-hand matrix, r x s.

 

K

target matrix, (p * r) x (q * s).

 

Returns

0 on success, E_NONCONF if matrix K is not correctly dimensioned for the operation.


gretl_matrix_kronecker_product_new ()

gretl_matrix *
gretl_matrix_kronecker_product_new (const gretl_matrix *A,
                                    const gretl_matrix *B,
                                    int *err);

Parameters

A

left-hand matrix, p x q.

 

B

right-hand matrix, r x s.

 

err

location to receive error code.

 

Returns

A newly allocated (p * r) x (q * s) matrix which is the Kronecker product of matrices A and B , or NULL on failure.


gretl_matrix_hdproduct ()

int
gretl_matrix_hdproduct (const gretl_matrix *A,
                        const gretl_matrix *B,
                        gretl_matrix *C);

Writes into C the horizontal direct product of A and B . That is, $C_i' = A_i' \otimes B_i'$ (in TeX notation). If B is NULL, then it's understood to be equal to A .

Parameters

A

left-hand matrix, r x p.

 

B

right-hand matrix, r x q or NULL.

 

C

target matrix, r x (p * q).

 

Returns

0 on success, E_NONCONF if A and B have different numbers of rows or matrix C is not correctly dimensioned for the operation.


gretl_matrix_hdproduct_new ()

gretl_matrix *
gretl_matrix_hdproduct_new (const gretl_matrix *A,
                            const gretl_matrix *B,
                            int *err);

If B is NULL, then it is implicitly taken as equal to A ; in this case, the returned matrix only contains the non-redundant elements; therefore, it has ncols = p*(p+1)/2 elements. Otherwise, all the products are computed and ncols = p*q.

Parameters

A

left-hand matrix, r x p.

 

B

right-hand matrix, r x q or NULL.

 

err

location to receive error code.

 

Returns

newly allocated r x ncols matrix which is the horizontal direct product of matrices A and B , or NULL on failure.


gretl_matrix_I_kronecker ()

int
gretl_matrix_I_kronecker (int p,
                          const gretl_matrix *B,
                          gretl_matrix *K);

Writes the Kronecker product of the identity matrix of order r and B into K .

Parameters

p

dimension of left-hand identity matrix.

 

B

right-hand matrix, r x s.

 

K

target matrix, (p * r) x (p * s).

 

Returns

0 on success, E_NONCONF if matrix K is not correctly dimensioned for the operation.


gretl_matrix_I_kronecker_new ()

gretl_matrix *
gretl_matrix_I_kronecker_new (int p,
                              const gretl_matrix *B,
                              int *err);

Writes the Kronecker product of the identity matrix of order r and B into a newly allocated matrix.

Parameters

p

dimension of left-hand identity matrix.

 

B

right-hand matrix, r x s.

 

err

location to receive error code.

 

Returns

the new matrix, or NULL on failure.


gretl_matrix_kronecker_I ()

int
gretl_matrix_kronecker_I (const gretl_matrix *A,
                          int r,
                          gretl_matrix *K);

Writes the Kronecker product of A and the identity matrix of order r into K .

Parameters

A

left-hand matrix, p x q.

 

r

dimension of right-hand identity matrix.

 

K

target matrix, (p * r) x (q * r).

 

Returns

0 on success, E_NONCONF if matrix K is not correctly dimensioned for the operation.


gretl_matrix_kronecker_I_new ()

gretl_matrix *
gretl_matrix_kronecker_I_new (const gretl_matrix *A,
                              int r,
                              int *err);

Writes into a newl allocated matrix the Kronecker product of A and the identity matrix of order r .

Parameters

A

left-hand matrix, p x q.

 

r

dimension of right-hand identity matrix.

 

err

location to receive error code.

 

Returns

the new matrix, or NULL on failure.


gretl_matrix_pow ()

gretl_matrix *
gretl_matrix_pow (const gretl_matrix *A,
                  double s,
                  int *err);

Calculates the matrix A^s. If s is a non-negative integer Golub and Van Loan's Algorithm 11.2.2 ("Binary Powering") is used. Otherwise A must be positive definite, and the power is computed via the eigen-decomposition of A .

Parameters

A

square source matrix.

 

s

exponent.

 

err

location to receive error code.

 

Returns

allocated matrix, or NULL on failure.


gretl_matrix_dot_product ()

double
gretl_matrix_dot_product (const gretl_matrix *a,
                          GretlMatrixMod amod,
                          const gretl_matrix *b,
                          GretlMatrixMod bmod,
                          int *err);

Parameters

a

left-hand matrix.

 

amod

modifier: GRETL_MOD_NONE or GRETL_MOD_TRANSPOSE.

 

b

right-hand matrix.

 

bmod

modifier: GRETL_MOD_NONE or GRETL_MOD_TRANSPOSE.

 

err

pointer to receive error code (zero on success, non-zero on failure), or NULL.

 

Returns

The dot (scalar) product of a (or a -transpose) and b (or b -transpose), or NADBL on failure.


gretl_vector_dot_product ()

double
gretl_vector_dot_product (const gretl_vector *a,
                          const gretl_vector *b,
                          int *err);

Parameters

a

first vector.

 

b

second vector.

 

err

pointer to receive error code (zero on success, non-zero on failure), or NULL.

 

Returns

The dot (scalar) product of a and b , or NADBL on failure.


gretl_rmatrix_vector_stat ()

gretl_matrix *
gretl_rmatrix_vector_stat (const gretl_matrix *m,
                           GretlVecStat vs,
                           int rowwise,
                           int skip_na,
                           int *err);

Parameters

m

source matrix.

 

vs

the required statistic or quantity: sum, product or mean.

 

rowwise

if non-zero go by rows, otherwise go by columns.

 

skip_na

ignore missing values.

 

err

location to receive error code.

 

Returns

a row vector or column vector containing the sums, products or means of the columns or rows of m . See also gretl_rmatrix_vector_stat() for the complex variant.


gretl_matrix_column_sd ()

gretl_matrix *
gretl_matrix_column_sd (const gretl_matrix *m,
                        int df,
                        int skip_na,
                        int *err);

Parameters

m

source matrix.

 

df

degrees of freedom for standard deviations.

 

skip_na

ignore missing values.

 

err

location to receive error code.

 

Returns

a row vector containing the standard deviations of the columns of m , or NULL on failure. If df is positive it is used as the divisor when calculating the column variance, otherwise the divisor is the number of rows in m .


gretl_matrix_demean_by_row ()

void
gretl_matrix_demean_by_row (gretl_matrix *m);

For each row of m , subtracts the row mean from each element on the row.

Parameters

m

matrix on which to operate.

 

gretl_matrix_standardize ()

int
gretl_matrix_standardize (gretl_matrix *m,
                          int dfcorr,
                          int skip_na);

Subtracts the column mean from each column of m and divides by the column standard deviation, using dfcorr as degrees of freedom correction (0 for MLE).

Parameters

m

matrix on which to operate.

 

dfcorr

degrees of freedom correction.

 

skip_na

ignore missing values.

 

gretl_matrix_center ()

int
gretl_matrix_center (gretl_matrix *m,
                     int skip_na);

Subtracts the column mean from each column of m .

Parameters

m

matrix on which to operate.

 

gretl_matrix_quantiles ()

gretl_matrix *
gretl_matrix_quantiles (const gretl_matrix *m,
                        const gretl_matrix *p,
                        int *err);

Parameters

m

matrix on which to operate.

 

p

vector of desired quantiles.

 

err

location to receive error code.

 

Returns

a matrix containing the p quantiles of the columns of m , or NULL on failure.


gretl_matrix_determinant ()

double
gretl_matrix_determinant (gretl_matrix *a,
                          int *err);

Compute the determinant of the square matrix a using the LU factorization. Matrix a is not preserved: it is overwritten by the factorization.

Parameters

a

gretl_matrix.

 

err

location to receive error code.

 

Returns

the determinant, or NADBL on failure.


gretl_matrix_log_determinant ()

double
gretl_matrix_log_determinant (gretl_matrix *a,
                              int *err);

Compute the log of the determinant of the square matrix a using LU factorization. Matrix a is not preserved: it is overwritten by the factorization.

Parameters

a

gretl_matrix.

 

err

location to receive error code.

 

Returns

the determinant, or NADBL on failure.


gretl_matrix_log_abs_determinant ()

double
gretl_matrix_log_abs_determinant (gretl_matrix *a,
                                  int *err);

Compute the log of the absolute value of the determinant of the square matrix a using LU factorization. Matrix a is not preserved: it is overwritten by the factorization.

Parameters

a

gretl_matrix.

 

err

location to receive error code.

 

Returns

the determinant, or NADBL on failure.


gretl_vcv_log_determinant ()

double
gretl_vcv_log_determinant (const gretl_matrix *m,
                           int *err);

Compute the log determinant of the symmetric positive-definite matrix m using Cholesky decomposition.

Parameters

m

gretl_matrix.

 

err

location to receive error code.

 

Returns

the log determinant, or NADBL on failure.


gretl_matrix_one_norm ()

double
gretl_matrix_one_norm (const gretl_matrix *m);

Parameters

m

gretl_matrix.

 

Returns

the 1-norm of m (the maximum value across the columns of m of the sum of the absolute values of the elements in the given column).


gretl_matrix_infinity_norm ()

double
gretl_matrix_infinity_norm (const gretl_matrix *m);

Parameters

m

gretl_matrix.

 

Returns

the infinity-norm of m (the maximum value across the rows of m of the sum of the absolute values of the elements in the given row).


gretl_LU_solve ()

int
gretl_LU_solve (gretl_matrix *a,
                gretl_matrix *b);

Solves ax = b for the unknown x, via LU decomposition using partial pivoting with row interchanges. On exit, b is replaced by the solution and a is replaced by its decomposition. Calls the LAPACK functions dgetrf() and dgetrs().

Parameters

a

square matrix.

 

b

matrix.

 

Returns

0 on successful completion, non-zero code on error.


gretl_LU_solve_invert ()

int
gretl_LU_solve_invert (gretl_matrix *a,
                       gretl_matrix *b);

Solves ax = b for the unknown x, using LU decomposition, then proceeds to use the decomposition to invert a . Calls the LAPACK functions dgetrf(), dgetrs() and dgetri(); the decomposition proceeds via partial pivoting with row interchanges.

On exit, b is replaced by the solution and a is replaced by its inverse.

Parameters

a

square matrix.

 

b

matrix.

 

Returns

0 on successful completion, non-zero code on error.


gretl_cholesky_decomp_solve ()

int
gretl_cholesky_decomp_solve (gretl_matrix *a,
                             gretl_matrix *b);

Solves ax = b for the unknown vector x, using Cholesky decomposition via the LAPACK functions dpotrf() and dpotrs().

On exit, b is replaced by the solution and a is replaced by its Cholesky decomposition.

Parameters

a

symmetric positive-definite matrix.

 

b

vector 'x' on input, solution 'b' on output.

 

Returns

0 on successful completion, or non-zero code on error.


gretl_cholesky_solve ()

int
gretl_cholesky_solve (const gretl_matrix *a,
                      gretl_vector *b);

Solves ax = b for the unknown vector x, using the pre-computed Cholesky decomposition of a . On exit, b is replaced by the solution.

Parameters

a

Cholesky-decomposed symmetric positive-definite matrix.

 

b

vector 'x'.

 

Returns

0 on successful completion, or non-zero code on error.


gretl_cholesky_invert ()

int
gretl_cholesky_invert (gretl_matrix *a);

Inverts a , which must contain the pre-computed Cholesky decomposition of a p.d. matrix, as may be obtained using gretl_cholesky_decomp_solve(). For speed there's no error checking of the input -- the caller should make sure it's OK.

Parameters

a

Cholesky-decomposed symmetric positive-definite matrix.

 

Returns

0 on successful completion, or non-zero code on error.


gretl_toeplitz_solve ()

gretl_vector *
gretl_toeplitz_solve (const gretl_vector *c,
                      const gretl_vector *r,
                      const gretl_vector *b,
                      double *det,
                      int *err);

Solves Tx = b for the unknown vector x, where T is a Toeplitz matrix, that is (zero-based)

T_{ij} = c_{i-j} for i <= j T_{ij} = r_{i-j} for i > j

Note that c[0] should equal r[0].

Parameters

c

Toeplitz column.

 

r

Toeplitz row.

 

b

right-hand side vector.

 

det

determinant (on exit)

 

err

error code.

 

Returns

a newly allocated vector, containing the solution, x.


gretl_matrix_XTX_new ()

gretl_matrix *
gretl_matrix_XTX_new (const gretl_matrix *X);

Parameters

X

matrix to process.

 

Returns

a newly allocated matrix containing X'X, or NULL on error.


gretl_inverse_from_cholesky_decomp ()

int
gretl_inverse_from_cholesky_decomp (gretl_matrix *targ,
                                    const gretl_matrix *src);

Computes in targ the inverse of a symmetric positive definite matrix, on the assumption that the original matrix this has already been Cholesky-decomposed in src .

Parameters

targ

matrix to hold inverse.

 

src

Cholesky-decomposed matrix.

 

Returns

0 on success; non-zero error code on failure.


gretl_invert_general_matrix ()

int
gretl_invert_general_matrix (gretl_matrix *a);

Computes the inverse of a general matrix using LU factorization. On exit a is overwritten with the inverse. Uses the LAPACK functions dgetrf() and dgetri().

Parameters

a

matrix to invert.

 

Returns

0 on success; non-zero error code on failure.


gretl_invert_symmetric_indef_matrix ()

int
gretl_invert_symmetric_indef_matrix (gretl_matrix *a);

Computes the inverse of a real symmetric matrix via the Bunch-Kaufman diagonal pivoting method. Uses the LAPACK functions dsytrf() and dsytri(). On exit a is overwritten with the inverse.

Parameters

a

matrix to invert.

 

Returns

0 on success; non-zero error code on failure.


gretl_invert_symmetric_matrix ()

int
gretl_invert_symmetric_matrix (gretl_matrix *a);

Computes the inverse of a symmetric positive definite matrix using Cholesky factorization. On exit a is overwritten with the inverse. Uses the LAPACK functions dpotrf() and dpotri().

Parameters

a

matrix to invert.

 

Returns

0 on success; non-zero error code on failure.


gretl_invert_symmetric_matrix2 ()

int
gretl_invert_symmetric_matrix2 (gretl_matrix *a,
                                double *ldet);

Computes the inverse of a symmetric positive definite matrix using Cholesky factorization, computing the log-determinant in the process. On exit a is overwritten with the inverse and if ldet is not NULL the log-determinant is written to that location. Uses the LAPACK functions dpotrf() and dpotri().

Parameters

a

matrix to invert.

 

ldet

location to receive log determinant, or NULL.

 

Returns

0 on success; non-zero error code on failure.


gretl_invert_packed_symmetric_matrix ()

int
gretl_invert_packed_symmetric_matrix (gretl_matrix *v);

Computes the inverse of a symmetric positive definite matrix, stored in vech form, using Cholesky factorization. On exit v is overwritten with the lower triangle of the inverse. Uses the LAPACK functions dpptrf() and dpptri().

Parameters

v

symmetric matrix in vech form (lower triangle packed as a column vector).

 

Returns

0 on success; non-zero error code on failure.


gretl_invert_triangular_matrix ()

int
gretl_invert_triangular_matrix (gretl_matrix *a,
                                char uplo);

Computes the inverse of a triangular matrix. On exit a is overwritten with the inverse. Uses the lapack function dtrtri.

Parameters

a

triangular matrix to invert.

 

uplo

'L' for lower triangular a , 'U' for upper.

 

Returns

0 on success; non-zero error code on failure.


gretl_invert_diagonal_matrix ()

int
gretl_invert_diagonal_matrix (gretl_matrix *a);

Computes the inverse of a diagonal matrix. On exit a is overwritten with the inverse.

Parameters

a

matrix to invert.

 

Returns

0 on success; non-zero error code on failure.


gretl_invert_matrix ()

int
gretl_invert_matrix (gretl_matrix *a);

Computes the inverse of matrix a : on exit a is overwritten with the inverse. If a is diagonal or symmetric, appropriate simple inversion routines are called.

Parameters

a

matrix to invert.

 

Returns

0 on success; non-zero error code on failure.


gretl_matrix_moore_penrose ()

int
gretl_matrix_moore_penrose (gretl_matrix *A,
                            double tol);

Computes the generalized inverse of matrix a via its SVD factorization, with the help of the lapack function dgesvd. On exit the original matrix is overwritten by the inverse.

Parameters

a

m x n matrix.

 

Returns

0 on success; non-zero error code on failure.


gretl_SVD_invert_matrix ()

int
gretl_SVD_invert_matrix (gretl_matrix *a);

Computes the inverse (or generalized inverse) of a general square matrix using SVD factorization, with the help of the lapack function dgesvd. If any of the singular values of a are less than 1.0e-9 the Moore-Penrose generalized inverse is computed instead of the standard inverse. On exit the original matrix is overwritten by the inverse.

Parameters

a

n x n matrix to invert.

 

Returns

0 on success; non-zero error code on failure.


gretl_invpd ()

int
gretl_invpd (gretl_matrix *a);

Computes the inverse of a symmetric positive definite matrix using Cholesky factorization. On exit a is overwritten with the inverse. Uses the LAPACK functions dpotrf() and dpotri(). Little checking is done, for speed: we assume the caller knows what he's doing.

Parameters

a

matrix to invert.

 

Returns

0 on success; non-zero error code on failure.


gretl_matrix_SVD ()

int
gretl_matrix_SVD (const gretl_matrix *a,
                  gretl_matrix **pu,
                  gretl_vector **ps,
                  gretl_matrix **pvt,
                  int full);

Computes SVD factorization of a general matrix using one of the the lapack functions dgesvd() or dgesdd(). A = U * diag(s) * Vt.

Parameters

x

m x n matrix to decompose.

 

pu

location for matrix U, or NULL if not wanted.

 

ps

location for vector of singular values, or NULL if not wanted.

 

pvt

location for matrix V (transposed), or NULL if not wanted.

 

full

if U and/or V are to be computed, a non-zero value flags production of "full-size" U (m x m) and/or V (n x n). Otherwise U will be m x min(m,n) and and V' will be min(m,n) x n. Note that this flag matters only if x is not square.

 

Returns

0 on success; non-zero error code on failure.


gretl_symmetric_matrix_rcond ()

double
gretl_symmetric_matrix_rcond (const gretl_matrix *m,
                              int *err);

Estimates the reciprocal condition number of the real symmetric positive definite matrix m (in the 1-norm), using the LAPACK functions dpotrf() and dpocon().

Parameters

m

matrix to examine.

 

err

location to receive error code.

 

Returns

the estimate, or NADBL on failure to allocate memory.


gretl_matrix_rcond ()

double
gretl_matrix_rcond (const gretl_matrix *m,
                    int *err);

Estimates the reciprocal condition number of the real matrix m (in the 1-norm).

Parameters

m

matrix to examine.

 

err

location to receive error code.

 

Returns

the estimate, or NADBL on failure to allocate memory.


gretl_matrix_cond_index ()

double
gretl_matrix_cond_index (const gretl_matrix *m,
                         int *err);

Estimates the condition number (a la Belsley) of the real matrix m .

Parameters

m

matrix to examine.

 

err

location to receive error code.

 

Returns

the estimate, or NADBL on failure.


gretl_symmetric_eigen_sort ()

int
gretl_symmetric_eigen_sort (gretl_matrix *evals,
                            gretl_matrix *evecs,
                            int rank);

Sorts the eigenvalues in evals from largest to smallest, and rearranges the columns in evecs correspondingly. If rank is greater than zero and less than the number of columns in evecs , then on output evecs is shrunk so that it contains only the columns associated with the largest rank eigenvalues.

Parameters

evals

array of real eigenvalues from symmetric matrix.

 

evecs

matrix of eigenvectors.

 

rank

desired number of columns in output.

 

Returns

0 on success; non-zero error code on failure.


gretl_general_matrix_eigenvals ()

gretl_matrix *
gretl_general_matrix_eigenvals (const gretl_matrix *m,
                                int *err);

Computes the eigenvalues of the general matrix m .

Parameters

m

square matrix on which to operate.

 

err

location to receive error code.

 

Returns

allocated matrix containing the eigenvalues, or NULL on failure. The returned matrix, on successful completion, is n x 2 (where n = the number of rows and columns in the matrix m ); the first column holds the real parts of the eigenvalues of m and the second the imaginary parts.


gretl_symmetric_matrix_eigenvals ()

gretl_matrix *
gretl_symmetric_matrix_eigenvals (gretl_matrix *m,
                                  int eigenvecs,
                                  int *err);

Computes the eigenvalues of the real symmetric matrix m . If eigenvecs is non-zero, also compute the orthonormal eigenvectors of m , which are stored in m . Uses the lapack function dsyevr(), or dsyev() for small matrices.

Parameters

m

n x n matrix to operate on.

 

eigenvecs

non-zero to calculate eigenvectors, 0 to omit.

 

err

location to receive error code.

 

Returns

n x 1 matrix containing the eigenvalues in ascending order, or NULL on failure.


gretl_symmetric_matrix_min_eigenvalue ()

double
gretl_symmetric_matrix_min_eigenvalue (const gretl_matrix *m);

gretl_symm_matrix_eigenvals_descending ()

gretl_matrix *
gretl_symm_matrix_eigenvals_descending
                               (gretl_matrix *m,
                                int eigenvecs,
                                int *err);

Computes the eigenvalues of the real symmetric matrix m . If eigenvecs is non-zero, also compute the orthonormal eigenvectors of m , which are stored in m . Uses the lapack function dsyev().

Parameters

m

n x n matrix to operate on.

 

eigenvecs

non-zero to calculate eigenvectors, 0 to omit.

 

err

location to receive error code.

 

Returns

n x 1 matrix containing the eigenvalues in descending order, or NULL on failure.


gretl_gensymm_eigenvals ()

gretl_matrix *
gretl_gensymm_eigenvals (const gretl_matrix *A,
                         const gretl_matrix *B,
                         gretl_matrix *V,
                         int *err);

Solves the generalized eigenvalue problem | A - \lambda B | = 0 , where both A and B are symmetric and B is positive definite.

Parameters

A

symmetric matrix.

 

B

symmetric positive definite matrix.

 

V

matrix to hold the generalized eigenvectors, or NULL if these are not required.

 

err

location to receive error code.

 

Returns

allocated storage containing the eigenvalues, in ascending order, or NULL on failure.


gretl_dgeev ()

gretl_matrix *
gretl_dgeev (const gretl_matrix *A,
             gretl_matrix *VR,
             gretl_matrix *VL,
             int *err);

old_eigengen ()

gretl_matrix *
old_eigengen (const gretl_matrix *m,
              gretl_matrix *VR,
              gretl_matrix *VL,
              int *err);

gretl_symm_matrix_lambda_min ()

double
gretl_symm_matrix_lambda_min (const gretl_matrix *m,
                              int *err);

Parameters

m

n x n matrix to operate on.

 

err

location to receive error code.

 

Returns

the minimum eigenvalue of the real symmetric matrix m , or NaN on error.


gretl_symm_matrix_lambda_max ()

double
gretl_symm_matrix_lambda_max (const gretl_matrix *m,
                              int *err);

Parameters

m

n x n matrix to operate on.

 

err

location to receive error code.

 

Returns

the maximum eigenvalue of the real symmetric matrix m , or NaN on error.


gretl_matrix_right_nullspace ()

gretl_matrix *
gretl_matrix_right_nullspace (const gretl_matrix *M,
                              int *err);

Given an m x n matrix M , construct a conformable matrix R such that MR = 0 (that is, all the columns of R are orthogonal to the space spanned by the rows of M ).

Parameters

M

matrix to operate on.

 

err

location to receive error code.

 

Returns

the allocated matrix R, or NULL on failure.


gretl_matrix_left_nullspace ()

gretl_matrix *
gretl_matrix_left_nullspace (const gretl_matrix *M,
                             GretlMatrixMod mod,
                             int *err);

Given an m x n matrix M , construct a conformable matrix L such that LM = 0 (that is, all the columns of M are orthogonal to the space spanned by the rows of L).

Parameters

M

matrix to operate on.

 

mod

GRETL_MOD_NONE or GRETL_MOD_TRANSPOSE

 

err

location to receive error code.

 

Returns

the allocated matrix L, or if mod is GRETL_MOD_TRANSPOSE, L', or NULL on failure.


gretl_matrix_row_concat ()

gretl_matrix *
gretl_matrix_row_concat (const gretl_matrix *a,
                         const gretl_matrix *b,
                         int *err);

Parameters

a

upper source matrix (m x n).

 

b

lower source matrix (p x n).

 

err

location to receive error code.

 

Returns

newly allocated matrix ((m+p) x n) that results from the row-wise concatenation of a and b , or NULL on failure.


gretl_matrix_col_concat ()

gretl_matrix *
gretl_matrix_col_concat (const gretl_matrix *a,
                         const gretl_matrix *b,
                         int *err);

Parameters

a

left-hand source matrix (m x n).

 

b

right-hand source matrix (m x p).

 

err

location to receive error code.

 

Returns

newly allocated matrix (m x (n+p)) that results from the column-wise concatenation of a and b , or NULL on failure.


gretl_matrix_direct_sum ()

gretl_matrix *
gretl_matrix_direct_sum (const gretl_matrix *a,
                         const gretl_matrix *b,
                         int *err);

Parameters

a

top left matrix.

 

b

bottom right matrix.

 

err

location to receive error code.

 

Returns

a new matrix containing the direct sum of a and b , or NULL on failure.


gretl_matrix_inplace_colcat ()

int
gretl_matrix_inplace_colcat (gretl_matrix *a,
                             const gretl_matrix *b,
                             const char *mask);

Concatenates onto a the selected columns of b , if the two matrices are conformable.

Parameters

a

matrix to be enlarged (m x n).

 

b

matrix from which columns should be added (m x p).

 

mask

char array, of length p, with 1s in positions corresponding to columns of b that are to be added to a , 0s elsewhere; or NULL to add all columns of b .

 

Returns

0 on success, non-zero code on error.


gretl_matrix_cumcol ()

gretl_matrix *
gretl_matrix_cumcol (const gretl_matrix *m,
                     int *err);

Parameters

m

source matrix.

 

err

error code.

 

Returns

a matrix of the same dimensions as m , containing the cumulated columns of m .


gretl_matrix_diffcol ()

gretl_matrix *
gretl_matrix_diffcol (const gretl_matrix *m,
                      double missval,
                      int *err);

Parameters

m

source matrix.

 

missval

value to represent missing observations.

 

err

error code.

 

Returns

a matrix of the same dimensions as m , containing missval in the first row and the difference between consecutive rows of m afterwards.


gretl_matrix_lag ()

gretl_matrix *
gretl_matrix_lag (const gretl_matrix *m,
                  const gretl_vector *k,
                  gretlopt opt,
                  double missval);

Parameters

m

source matrix.

 

k

vector of lag orders (> 0 for lags, < 0 for leads).

 

opt

use OPT_L to arrange multiple lags by lag rather than by variable.

 

missval

value to represent missing observations.

 

Returns

A matrix of the same dimensions as m , containing lags of the variables in the columns of m , with missing values set to missval .


gretl_matrix_inplace_lag ()

int
gretl_matrix_inplace_lag (gretl_matrix *targ,
                          const gretl_matrix *src,
                          int k);

Fills out targ (if it is of the correct dimensions), with (columnwise) lags of src , using 0 for missing values.

Parameters

targ

target matrix.

 

src

source matrix.

 

k

lag order (> 0 for lags, < 0 for leads).

 

Returns

0 on success, non-zero code otherwise.


gretl_matrix_cholesky_decomp ()

int
gretl_matrix_cholesky_decomp (gretl_matrix *a);

Computes the Cholesky factorization of the symmetric, positive definite matrix a . On exit the lower triangle of a is replaced by the factor L, as in a = LL', and the upper triangle is set to zero. Uses the lapack function dpotrf.

Parameters

a

matrix to operate on.

 

Returns

0 on success; 1 on failure.


gretl_matrix_psd_root ()

int
gretl_matrix_psd_root (gretl_matrix *a,
                       int check);

Computes the LL' factorization of the symmetric, positive semidefinite matrix a . On successful exit the lower triangle of a is replaced by the factor L and the upper triangle is set to zero.

Parameters

a

matrix to operate on.

 

check

if non-zero, perform a test for psd status.

 

Returns

0 on success; non-zero on failure.


gretl_matrix_QR_decomp ()

int
gretl_matrix_QR_decomp (gretl_matrix *M,
                        gretl_matrix *R);

Computes the QR factorization of M . On successful exit the matrix M holds Q, and, if R is not NULL, the upper triangle of R holds R. Uses the LAPACK functions dgeqrf() and dorgqr().

Parameters

M

m x n matrix to be decomposed.

 

R

n x n matrix into which to write R, as in M = Q * R, or NULL if this is not wanted.

 

Returns

0 on success, non-zero on failure.


gretl_matrix_QR_pivot_decomp ()

int
gretl_matrix_QR_pivot_decomp (gretl_matrix *M,
                              gretl_matrix *R,
                              gretl_matrix *P);

gretl_triangular_matrix_rcond ()

double
gretl_triangular_matrix_rcond (const gretl_matrix *A,
                               char uplo,
                               char diag);

Parameters

A

triangular matrix.

 

uplo

'U' for upper triangular input, 'L' for lower.

 

diag

'N' for non-unit triangular, 'U' for unit triangular.

 

Returns

the reciprocal condition number of R in the 1-norm, or NADBL on failure.


gretl_check_QR_rank ()

int
gretl_check_QR_rank (const gretl_matrix *R,
                     int *err,
                     double *rcnd);

Checks the reciprocal condition number of R and calculates the rank of the matrix QR. If rcnd is not NULL it receives the reciprocal condition number.

Parameters

R

matrix R from QR decomposition.

 

err

location to receive error code.

 

rcnd

location to receive reciprocal condition number.

 

Returns

on success, the rank of QR.


gretl_matrix_rank ()

int
gretl_matrix_rank (const gretl_matrix *a,
                   double eps,
                   int *err);

Computes the rank of a via its SV decomposition.

Parameters

a

matrix to examine.

 

eps

value below which a singular value is taken to be zero. Give 0 or NADBL for automatic use of machine epsilon.

 

err

location to receive error code on failure.

 

Returns

the rank of a , or 0 on failure.


gretl_matrix_ols ()

int
gretl_matrix_ols (const gretl_vector *y,
                  const gretl_matrix *X,
                  gretl_vector *b,
                  gretl_matrix *vcv,
                  gretl_vector *uhat,
                  double *s2);

Computes OLS estimates using Cholesky factorization by default, but with a fallback to QR decomposition if the data are highly ill-conditioned, and puts the coefficient estimates in b . Optionally, calculates the covariance matrix in vcv and the residuals in uhat .

Parameters

y

dependent variable vector.

 

X

matrix of independent variables.

 

b

vector to hold coefficient estimates.

 

vcv

matrix to hold the covariance matrix of the coefficients, or NULL if this is not needed.

 

uhat

vector to hold the regression residuals, or NULL if these are not needed.

 

s2

pointer to receive residual variance, or NULL. Note: if s2 is NULL, the "vcv" estimate will be plain (X'X)^{-1}.

 

Returns

0 on success, non-zero error code on failure.


gretl_matrix_multi_ols ()

int
gretl_matrix_multi_ols (const gretl_matrix *Y,
                        const gretl_matrix *X,
                        gretl_matrix *B,
                        gretl_matrix *E,
                        gretl_matrix **XTXi);

Computes OLS estimates using Cholesky decomposition by default, but with a fallback to QR decomposition if the data are highly ill-conditioned, and puts the coefficient estimates in B . Optionally, calculates the residuals in E .

Parameters

y

T x g matrix of dependent variables.

 

X

T x k matrix of independent variables.

 

B

k x g matrix to hold coefficient estimates, or NULL.

 

E

T x g matrix to hold the regression residuals, or NULL if these are not needed.

 

XTXi

location to receive (X'X)^{-1} on output, or NULL if this is not needed.

 

Returns

0 on success, non-zero error code on failure.


gretl_matrix_multi_SVD_ols ()

int
gretl_matrix_multi_SVD_ols (const gretl_matrix *Y,
                            const gretl_matrix *X,
                            gretl_matrix *B,
                            gretl_matrix *E,
                            gretl_matrix **XTXi);

Computes OLS estimates using SVD decomposition, and puts the coefficient estimates in B . Optionally, calculates the residuals in E , (X'X)^{-1} in XTXi .

Parameters

y

T x g matrix of dependent variables.

 

X

T x k matrix of independent variables.

 

B

k x g matrix to hold coefficient estimates, or NULL.

 

E

T x g matrix to hold the regression residuals, or NULL if these are not needed.

 

XTXi

location to receive (X'X)^{-1}, or NULL if this is not needed.

 

Returns

0 on success, non-zero error code on failure.


gretl_matrix_QR_ols ()

int
gretl_matrix_QR_ols (const gretl_matrix *Y,
                     const gretl_matrix *X,
                     gretl_matrix *B,
                     gretl_matrix *E,
                     gretl_matrix **XTXi,
                     gretl_matrix **Qout);

Computes OLS estimates using QR decomposition, and puts the coefficient estimates in B . Optionally, calculates the residuals in E , (X'X)^{-1} in XTXi , and/or the matrix Q in Qout .

Parameters

y

T x g matrix of dependent variables.

 

X

T x k matrix of independent variables.

 

B

k x g matrix to hold coefficient estimates.

 

E

T x g matrix to hold the regression residuals, or NULL if these are not needed.

 

XTXi

location to receive (X'X)^{-1}, or NULL if this is not needed.

 

Qout

location to receive Q on output, or NULL.

 

Returns

0 on success, non-zero error code on failure.


gretl_matrix_r_squared ()

double
gretl_matrix_r_squared (const gretl_matrix *y,
                        const gretl_matrix *X,
                        const gretl_matrix *b,
                        int *err);

Parameters

y

dependent variable, T-vector.

 

X

independent variables matrix, T x k.

 

b

coefficients, k-vector.

 

err

location to receive error code.

 

Returns

the unadjusted R-squared, based on the regression represented by y , X and b , or NADBL on failure.


gretl_matrix_SVD_johansen_solve ()

int
gretl_matrix_SVD_johansen_solve (const gretl_matrix *R0,
                                 const gretl_matrix *R1,
                                 gretl_matrix *evals,
                                 gretl_matrix *B,
                                 gretl_matrix *A,
                                 int jrank);

Solves the Johansen generalized eigenvalue problem via SVD decomposition. See J. A. Doornik and R. J. O'Brien, "Numerically stable cointegration analysis", Computational Statistics and Data Analysis, 41 (2002), pp. 185-193, Algorithm 4.

If B is non-null it should be p1 x p on input; it will be trimmed to p1 x jrank on output if jrank < p. If A is non-null it should be p x p on input; it will be trimmed to p x jrank on output if jrank < p. evals should be a vector of length jrank .

Parameters

R0

T x p matrix of residuals.

 

R1

T x p1 matrix of residuals.

 

evals

vector to receive eigenvals, or NULL if not wanted.

 

B

matrix to hold \beta, or NULL if not wanted.

 

A

matrix to hold \alpha, or NULL if not wanted.

 

jrank

cointegration rank, <= p.

 

Returns

0 on success; non-zero error code on failure.


gretl_matrix_restricted_ols ()

int
gretl_matrix_restricted_ols (const gretl_vector *y,
                             const gretl_matrix *X,
                             const gretl_matrix *R,
                             const gretl_vector *q,
                             gretl_vector *b,
                             gretl_matrix *vcv,
                             gretl_vector *uhat,
                             double *s2);

Computes OLS estimates restricted by R and q, using the lapack function dgglse(), and puts the coefficient estimates in b . Optionally, calculates the covariance matrix in vcv . If q is NULL this is taken as equivalent to a zero vector.

Parameters

y

dependent variable vector.

 

X

matrix of independent variables.

 

R

left-hand restriction matrix, as in Rb = q.

 

q

right-hand restriction vector or NULL.

 

b

vector to hold coefficient estimates.

 

vcv

matrix to hold the covariance matrix of the coefficients, or NULL if this is not needed.

 

uhat

vector to hold residuals, if wanted.

 

s2

pointer to receive residual variance, or NULL. If vcv is non-NULL and s2 is NULL, the vcv estimate is just W^{-1}.

 

Returns

0 on success, non-zero error code on failure.


gretl_matrix_restricted_multi_ols ()

int
gretl_matrix_restricted_multi_ols (const gretl_matrix *Y,
                                   const gretl_matrix *X,
                                   const gretl_matrix *R,
                                   const gretl_matrix *q,
                                   gretl_matrix *B,
                                   gretl_matrix *U,
                                   gretl_matrix **W);

Computes LS estimates restricted by R and q , putting the coefficient estimates into B . The R matrix must have g*k columns; each row represents a linear restriction; q must be a column vector of length equal to the number of rows of R .

Parameters

Y

dependent variable matrix, T x g.

 

X

matrix of independent variables, T x k.

 

R

left-hand restriction matrix, as in RB = q.

 

q

right-hand restriction matrix.

 

B

matrix to hold coefficient estimates, k x g.

 

U

matrix to hold residuals (T x g), if wanted.

 

pW

pointer to matrix to hold the RLS counterpart to (X'X)^{-1}, if wanted.

 

Returns

0 on success, non-zero error code on failure.


gretl_matrix_SVD_ols ()

int
gretl_matrix_SVD_ols (const gretl_vector *y,
                      const gretl_matrix *X,
                      gretl_vector *b,
                      gretl_matrix *vcv,
                      gretl_vector *uhat,
                      double *s2);

Computes OLS estimates using SVD decomposition, and puts the coefficient estimates in b . Optionally, calculates the covariance matrix in vcv and the residuals in uhat .

Parameters

y

dependent variable vector.

 

X

matrix of independent variables.

 

b

vector to hold coefficient estimates.

 

vcv

matrix to hold the covariance matrix of the coefficients, or NULL if this is not needed.

 

uhat

vector to hold the regression residuals, or NULL if these are not needed.

 

s2

pointer to receive residual variance, or NULL. Note: if s2 is NULL, the vcv estimate will be plain (X'X)^{-1}.

 

Returns

0 on success, non-zero error code on failure.


gretl_matrix_qform ()

int
gretl_matrix_qform (const gretl_matrix *A,
                    GretlMatrixMod amod,
                    const gretl_matrix *X,
                    gretl_matrix *C,
                    GretlMatrixMod cmod);

Computes either A * X * A' (if amod = GRETL_MOD_NONE) or A' * X * A (if amod = GRETL_MOD_TRANSPOSE), with the result written into C . The matrix X must be symmetric, but this is not checked, to save time. If you are in doubt on this point you can call gretl_matrix_is_symmetric() first.

Parameters

A

m * k matrix or k * m matrix, depending on amod .

 

amod

GRETL_MOD_NONE or GRETL_MOD_TRANSPOSE: in the first case A should be m * k; in the second, k * m;

 

X

symmetric k * k matrix.

 

C

matrix to hold the product.

 

cmod

modifier: GRETL_MOD_NONE, or GRETL_MOD_CUMULATE to add the result to the existing value of C , or GRETL_MOD_DECREMENT to subtract from the existing value of C .

 

Returns

0 on success; non-zero error code on failure.


gretl_matrix_diag_qform ()

int
gretl_matrix_diag_qform (const gretl_matrix *A,
                         GretlMatrixMod amod,
                         const gretl_vector *X,
                         gretl_matrix *C,
                         GretlMatrixMod cmod);

Computes either A * md * A' (if amod = GRETL_MOD_NONE) or A' * md * A (if amod = GRETL_MOD_TRANSPOSE), where md is a diagonal matrix holding the vector d . The result is written into C .

Parameters

A

m * k matrix or k * m matrix, depending on amod .

 

amod

GRETL_MOD_NONE or GRETL_MOD_TRANSPOSE: in the first case A should be m * k; in the second, k * m;

 

d

k-element vector.

 

C

matrix to hold the product.

 

cmod

modifier: GRETL_MOD_NONE, or GRETL_MOD_CUMULATE to add the result to the existing value of C , or GRETL_MOD_DECREMENT to subtract from the existing value of C .

 

Returns

0 on success; non-zero error code on failure.


gretl_scalar_qform ()

double
gretl_scalar_qform (const gretl_vector *b,
                    const gretl_matrix *X,
                    int *err);

Computes the scalar product bXb', or b'Xb if b is a column vector. The content of err is set to a non-zero code on failure.

Parameters

b

k-vector.

 

X

symmetric k x k matrix.

 

err

pointer to receive error code.

 

Returns

the scalar product, or NADBL on failure.


gretl_matrix_columnwise_product ()

int
gretl_matrix_columnwise_product (const gretl_matrix *A,
                                 const gretl_matrix *B,
                                 const gretl_matrix *S,
                                 gretl_matrix *C);

gretl_matrix_diagonal_sandwich ()

int
gretl_matrix_diagonal_sandwich (const gretl_vector *d,
                                const gretl_matrix *X,
                                gretl_matrix *DXD);

Computes in DXD (which must be pre-allocated), the matrix product D * X * D, where D is a diagonal matrix with elements given by the vector d .

Parameters

d

k-vector.

 

X

k * k matrix.

 

DXD

target k * k matrix.

 

Returns

0 on success, non-zero code on error.


gretl_matrix_set_t1 ()

int
gretl_matrix_set_t1 (gretl_matrix *m,
                     int t);

Sets an integer value on m , which can be retrieved using gretl_matrix_get_t1().

Parameters

m

matrix to operate on.

 

t

integer value to set.

 

Returns

0 on success, non-ero on error.


gretl_matrix_set_t2 ()

int
gretl_matrix_set_t2 (gretl_matrix *m,
                     int t);

Sets an integer value on m , which can be retrieved using gretl_matrix_get_t2().

Parameters

m

matrix to operate on.

 

t

integer value to set.

 

Returns

0 on success, non-ero on error.


gretl_matrix_get_t1 ()

int
gretl_matrix_get_t1 (const gretl_matrix *m);

Parameters

m

matrix to read from.

 

Returns

the integer that has been set on m using gretl_matrix_set_t1(), or zero if no such value has been set.


gretl_matrix_get_t2 ()

int
gretl_matrix_get_t2 (const gretl_matrix *m);

Parameters

m

matrix to read from.

 

Returns

the integer that has been set on m using gretl_matrix_set_t2(), or zero if no such value has been set.


gretl_matrix_is_dated ()

int
gretl_matrix_is_dated (const gretl_matrix *m);

Parameters

m

matrix to examine.

 

Returns

1 if matrix m has integer indices recorded via gretl_matrix_set_t1() and gretl_matrix_set_t2(), such that t1 >= 0 and t2 > t1, otherwise zero.


gretl_is_identity_matrix ()

int
gretl_is_identity_matrix (const gretl_matrix *m);

Parameters

m

matrix to examine.

 

Returns

1 if m is an identity matrix, 0 otherwise.


gretl_is_zero_matrix ()

int
gretl_is_zero_matrix (const gretl_matrix *m);

Parameters

m

matrix to examine.

 

Returns

1 if m is a zero matrix, 0 otherwise.


gretl_matrix_isfinite ()

gretl_matrix *
gretl_matrix_isfinite (const gretl_matrix *m,
                       int *err);

Parameters

m

matrix to examine.

 

err

location to receive error code.

 

Returns

a matrix with 1s in positions corresponding to finite elements of m , zeros otherwise.


gretl_matrix_get_structure ()

int
gretl_matrix_get_structure (const gretl_matrix *m);

gretl_matrices_are_equal ()

int
gretl_matrices_are_equal (const gretl_matrix *a,
                          const gretl_matrix *b,
                          double tol,
                          int *err);

Parameters

a

first matrix in comparison.

 

b

second matrix in comparison.

 

tol

numerical tolerance.

 

err

location to receive error code.

 

Returns

1 if the matrices a and b compare equal, 0 if they differ, and -1 if the comparison is invalid, in which case E_NONCONF is written to err .


gretl_covariance_matrix ()

gretl_matrix *
gretl_covariance_matrix (const gretl_matrix *m,
                         int corr,
                         int dfc,
                         int *err);

Parameters

m

(x x n) matrix containing n observations on each of k variables.

 

corr

flag for computing correlations.

 

dfc

degrees of freedom correction: use 1 for sample variance, 0 for MLE.

 

err

pointer to receive non-zero error code in case of failure, or NULL.

 

Returns

the covariance matrix of variables in the columns of m , or the correlation matrix if corr is non-zero.


gretl_matrix_GG_inverse ()

gretl_matrix *
gretl_matrix_GG_inverse (const gretl_matrix *G,
                         int *err);

Multiples G' into G and inverts the result. A shortcut function intended for producing an approximation to the Hessian given a gradient matrix.

Parameters

G

T x k source matrix.

 

err

location to receive error code.

 

Returns

the newly allocated k x k inverse on success, or NULL on error.


gretl_matrix_varsimul ()

gretl_matrix *
gretl_matrix_varsimul (const gretl_matrix *A,
                       const gretl_matrix *U,
                       const gretl_matrix *x0,
                       int *err);

Simulates a p-order n-variable VAR: x_t = \sum A_i x_{t-i} + u_t

The A_i matrices must be stacked horizontally into the A argument, that is: A = A_1 ~ A_2 ~ A_p. The u_t vectors are contained (as rows) in U . Initial values are in x0 .

Note the that the arrangement of the A matrix is somewhat sub-optimal computationally, since its elements have to be reordered by the function reorder_A (see above). However, the present form is more intuitive for a human being, and that's what counts.

Parameters

A

n x np coefficient matrix.

 

U

T x n data matrix.

 

x0

p x n matrix for initialization.

 

err

location to receive error code.

 

Returns

a newly allocated T+p x n matrix on success, whose t-th row is (x_t)', or NULL on error.


gretl_matrix_array_new ()

gretl_matrix **
gretl_matrix_array_new (int n);

Allocates an array of n gretl matrix pointers. On successful allocation of the array, each element is initialized to NULL.

Parameters

n

number of matrices.

 

Returns

pointer on sucess, NULL on failure.


gretl_matrix_array_new_with_size ()

gretl_matrix **
gretl_matrix_array_new_with_size (int n,
                                  int rows,
                                  int cols);

Allocates an array of n gretl matrix pointers, each one with size rows * cols .

Parameters

n

number of matrices.

 

rows

number of rows in each matrix.

 

cols

number of columns in each matrix.

 

Returns

pointer on sucess, NULL on failure.


gretl_matrix_array_free ()

void
gretl_matrix_array_free (gretl_matrix **A,
                         int n);

Frees each of the n gretl matrices in the array A , and the array itself. See also gretl_matrix_array_alloc().

Parameters

A

dyamically allocated array of gretl matrices.

 

n

number of matrices in array.

 

gretl_matrix_values ()

gretl_matrix *
gretl_matrix_values (const double *x,
                     int n,
                     gretlopt opt,
                     int *err);

Parameters

x

array to process.

 

n

length of array.

 

opt

if OPT_S the array of values will be sorted, otherwise given in order of occurrence.

 

err

location to receive error code.

 

Returns

an allocated matrix containing the distinct values in array x , or NULL on failure.


gretl_matrix_shape ()

gretl_matrix *
gretl_matrix_shape (const gretl_matrix *A,
                    int r,
                    int c,
                    int *err);

Creates an (r x c) matrix containing the re-arranged values of A. Elements are read from A by column and written into the target, also by column. If A contains less elements than n = r*c, they are repeated cyclically; if A has more elements, only the first n are used.

Parameters

A

array to process.

 

r

rows of target matrix.

 

c

columns of target matrix.

 

Returns

the generated matrix, or NULL on failure.


gretl_matrix_trim_rows ()

gretl_matrix *
gretl_matrix_trim_rows (const gretl_matrix *A,
                        int ttop,
                        int tbot,
                        int *err);

Creates a new matrix which is a copy of A with ttop rows trimmed from the top and tbot rows trimmed from the bottom.

Parameters

A

array to process.

 

ttop

rows to trim at top.

 

tbot

rows to trim at bottom.

 

err

location to receive error code.

 

Returns

the generated matrix, or NULL on failure.


gretl_matrix_minmax ()

gretl_matrix *
gretl_matrix_minmax (const gretl_matrix *A,
                     int mm,
                     int rc,
                     int idx,
                     int skip_na,
                     int *err);

Creates a matrix holding the row or column mimima or maxima from A , either as values or as location indices. For example, if mm = 0, rc = 0, and idx = 0, the created matrix is m x 1 and holds the values of the row minima.

Parameters

A

m x n matrix to examine.

 

mm

0 for minima, 1 for maxima.

 

rc

0 for row-wise, 1 for column-wise.

 

idx

0 for values, 1 for indices.

 

skip_na

ignore missing values.

 

err

location to receive error code.

 

Returns

the generated matrix, or NULL on failure.


gretl_matrix_global_minmax ()

double
gretl_matrix_global_minmax (const gretl_matrix *A,
                            int mm,
                            int *err);

Parameters

A

matrix to examine.

 

mm

0 for minimum, 1 for maximum.

 

err

location to receive error code.

 

Returns

the smallest or greatest element of A , ignoring NaNs but not infinities (?), or NADBL on failure.


gretl_matrix_global_sum ()

double
gretl_matrix_global_sum (const gretl_matrix *A,
                         int *err);

Parameters

A

matrix to examine.

 

err

location to receive error code.

 

Returns

the sum of the elements of A , or NADBL on failure.


gretl_matrix_pca ()

gretl_matrix *
gretl_matrix_pca (const gretl_matrix *X,
                  int p,
                  gretlopt opt,
                  int *err);

Carries out a Principal Components analysis of X and returns the first p components: the component corresponding to the largest eigenvalue of the correlation matrix of X is placed in column 1, and so on.

Parameters

X

T x m data matrix.

 

p

number of principal components to return: 0 < p <= m.

 

opt

if OPT_V, use the covariance matrix rather than the correlation matrix as basis.

 

err

location to receive error code.

 

Returns

the generated matrix, or NULL on failure.


gretl_matrix_xtab ()

gretl_matrix *
gretl_matrix_xtab (int t1,
                   int t2,
                   const double *x,
                   const double *y,
                   int *err);

Computes the cross tabulation of the values contained in the vectors x (by row) and y (by column). These must be integer values.

Parameters

x

data vector

 

y

data vector

 

t1

start

 

t2

end

 

err

error code

 

Returns

the generated matrix, or NULL on failure.


matrix_matrix_xtab ()

gretl_matrix *
matrix_matrix_xtab (const gretl_matrix *x,
                    const gretl_matrix *y,
                    int *err);

Computes the cross tabulation of the values contained in the vectors x (by row) and y (by column). These must be integer values.

Parameters

x

data vector

 

y

data vector

 

err

error code

 

Returns

the generated matrix, or NULL on failure.


gretl_matrix_bool_sel ()

gretl_matrix *
gretl_matrix_bool_sel (const gretl_matrix *A,
                       const gretl_matrix *sel,
                       int rowsel,
                       int *err);

If rowsel = 1, constructs a matrix which contains the rows of A corresponding to non-zero values in the vector sel ; if rowsel = 0, does the same thing but column-wise.

Parameters

A

matrix.

 

sel

selection vector.

 

rowsel

row/column mode selector.

 

err

location to receive error code.

 

Returns

the generated matrix, or NULL on failure.


gretl_matrix_sort_by_column ()

gretl_matrix *
gretl_matrix_sort_by_column (const gretl_matrix *m,
                             int k,
                             int *err);

Produces a matrix which contains the rows of m , re- ordered by increasing value of the elements in column k .

Parameters

m

matrix.

 

k

column by which to sort.

 

err

location to receive error code.

 

Returns

the generated matrix, or NULL on failure.


gretl_vector_sort ()

gretl_matrix *
gretl_vector_sort (const gretl_matrix *v,
                   int descending,
                   int *err);

gretl_matrix_covariogram ()

gretl_matrix *
gretl_matrix_covariogram (const gretl_matrix *X,
                          const gretl_matrix *u,
                          const gretl_matrix *w,
                          int p,
                          int *err);

Produces the matrix covariogram,

\sum_{j=-p}^{p} \sum_j w_{|j|} (X_t' u_t u_{t-j} X_{t-j})

If u is not given the u terms are omitted, and if w is not given, all the weights are 1.0.

Parameters

X

T x k matrix (typically containing regressors).

 

u

T-vector (typically containing residuals), or NULL.

 

w

(p+1)-vector of weights, or NULL.

 

p

lag order >= 0.

 

err

location to receive error code.

 

Returns

the generated matrix, or NULL on failure.


gretl_matrix_commute ()

gretl_matrix *
gretl_matrix_commute (gretl_matrix *A,
                      int r,
                      int c,
                      int pre,
                      int add_id,
                      int *err);

It is assumed that A is a matrix with (r *c ) rows, so each of its columns can be seen as the vectorization of an (r x c ) matrix. Each column of the output matrix contains the vectorization of the transpose of the corresponding column of A . This is equivalent to premultiplying A by the so-called "commutation matrix" $K_{r,c}$. If the add_id flag is non-zero, then A is added to the output matrix, so that A is premultiplied by (I + K_{r,c}) if pre is nonzero, postmultiplied if pre is 0.

See eg Magnus and Neudecker (1988), "Matrix Differential Calculus with Applications in Statistics and Econometrics".

Parameters

A

source matrix.

 

r

row dimension.

 

c

column dimension.

 

pre

premultiply (Boolean flag).

 

add_id

add identity matrix (Boolean flag).

 

err

location to receive error code.

 

gretl_matrix_transcribe_obs_info ()

void
gretl_matrix_transcribe_obs_info (gretl_matrix *targ,
                                  const gretl_matrix *src);

If targ and src have the same number of rows, and if the rows of src are identified by observation stamps while those of targ are not so identified, copy the stamp information across to targ . (Or if the given conditions are not satified, do nothing.)

Parameters

targ

target matrix.

 

src

source matrix.

 

gretl_matrix_set_colnames ()

int
gretl_matrix_set_colnames (gretl_matrix *m,
                           char **S);

Sets an array of strings on m which can be retrieved using gretl_matrix_get_colnames(). Note that S must contain as many strings as m has columns. The matrix takes ownership of S , which should be allocated and not subsequently touched by the caller.

Parameters

m

target matrix.

 

S

array of strings.

 

Returns

0 on success, non-zero code on error.


gretl_matrix_set_rownames ()

int
gretl_matrix_set_rownames (gretl_matrix *m,
                           char **S);

Sets an array of strings on m which can be retrieved using gretl_matrix_get_rownames(). Note that S must contain as many strings as m has rows. The matrix takes ownership of S , which should be allocated and not subsequently touched by the caller.

Parameters

m

target matrix.

 

S

array of strings.

 

Returns

0 on success, non-zero code on error.


gretl_matrix_get_colnames ()

const char **
gretl_matrix_get_colnames (const gretl_matrix *m);

Parameters

m

matrix

 

Returns

The array of strings set on m using gretl_matrix_set_colnames(), or NULL if no such strings have been set. The returned array will contain as many strings as m has columns.


gretl_matrix_get_rownames ()

const char **
gretl_matrix_get_rownames (const gretl_matrix *m);

Parameters

m

matrix

 

Returns

The array of strings set on m using gretl_matrix_set_rownames(), or NULL if no such strings have been set. The returned array will contain as many strings as m has rows.


gretl_matrix_destroy_info ()

void
gretl_matrix_destroy_info (gretl_matrix *m);

lapack_mem_free ()

void
lapack_mem_free (void);

Cleanup function, called by libgretl_cleanup(). Frees any memory that has been allocated internally as temporary workspace for LAPACK functions.


set_blas_mnk_min ()

void
set_blas_mnk_min (int mnk);

By default all matrix multiplication within libgretl is done using native code, but there is the possibility of farming out multiplication to the BLAS.

When multiplying an m x n matrix into an n x k matrix libgretl finds the product of the dimensions, m*n*k, and compares this with an internal threshhold variable, blas_mnk_min. If and only if blas_mnk_min >= 0 and n*m*k >= blas_mnk_min, then we use the BLAS. By default blas_mnk_min is set to -1 (BLAS never used).

If you have an optimized version of the BLAS you may want to set blas_mnk_min to some suitable positive value. (Setting it to 0 would result in external calls to the BLAS for all matrix multiplications, however small, which is unlikely to be optimal.)

Parameters

mnk

value to set.

 

get_blas_mnk_min ()

int
get_blas_mnk_min (void);

Returns

the value of the internal variable blas_mnk_min. See set_blas_mnk_min().


set_simd_k_max ()

void
set_simd_k_max (int k);

get_simd_k_max ()

int
get_simd_k_max (void);

set_simd_mn_min ()

void
set_simd_mn_min (int mn);

get_simd_mn_min ()

int
get_simd_mn_min (void);

Types and Values

R_DIAG_MIN

#define R_DIAG_MIN 1.0e-8

enum GretlMatrixMod

Members

GRETL_MOD_NONE

   

GRETL_MOD_TRANSPOSE

   

GRETL_MOD_SQUARE

   

GRETL_MOD_CUMULATE

   

GRETL_MOD_DECREMENT

   

GRETL_MOD_CTRANSP

   

enum GretlMatrixStructure

Members

GRETL_MATRIX_SQUARE

   

GRETL_MATRIX_LOWER_TRIANGULAR

   

GRETL_MATRIX_UPPER_TRIANGULAR

   

GRETL_MATRIX_SYMMETRIC

   

GRETL_MATRIX_DIAGONAL

   

GRETL_MATRIX_IDENTITY

   

GRETL_MATRIX_SCALAR

   

enum GretlVecStat

Members

V_SUM

   

V_PROD

   

V_MEAN

   

enum ConfType

Members

CONF_NONE

   

CONF_ELEMENTS

   

CONF_A_COLVEC

   

CONF_B_COLVEC

   

CONF_A_ROWVEC

   

CONF_B_ROWVEC

   

CONF_A_SCALAR

   

CONF_B_SCALAR

   

CONF_AC_BR

   

CONF_AR_BC

   

gretl_vector

typedef struct gretl_matrix_ gretl_vector;

matrix_info

typedef struct matrix_info_ matrix_info;

gretl_matrix

typedef struct {
    int rows;
    int cols;
    double *val;
    double _Complex *z; /* was "complex" */
    int is_complex;
} gretl_matrix;

The basic libgretl matrix type; gretl_vector is an alias that can be used for matrices with rows or cols = 1.

Members

int rows;

number of rows in matrix

 

int cols;

number of columns

 

double *val;

flat array of double-precision values

 

int is_complex;

   

gretl_matrix_block

typedef struct gretl_matrix_block_ gretl_matrix_block;