| Libgretl Reference Manual | ||||
|---|---|---|---|---|
| Top | Description | ||||
#define M_NA #define R_DIAG_MIN enum GretlMatrixMod; enum GretlMatrixStructure; gretl_matrix; typedef gretl_vector; gretl_matrix_block; #define gretl_matrix_get (m,i,j) #define gretl_vector_get (v,i) #define gretl_matrix_set (m,i,j,x) #define gretl_vector_set (v,i,x) #define gretl_matrix_cols (m) #define gretl_matrix_rows (m) #define gretl_vector_get_length (v) #define gretl_vector_alloc (i) #define gretl_column_vector_alloc (i) #define gretl_vector_free (v) #define gretl_matrix_is_scalar (m) #define gretl_is_null_matrix (m) int get_gretl_matrix_err (void); void clear_gretl_matrix_err (void); void gretl_matrix_print (const gretl_matrix *m, const char *msg); int gretl_matrix_xna_check (const gretl_matrix *m); int gretl_matrix_is_symmetric (const gretl_matrix *m); int gretl_matrix_is_idempotent (const gretl_matrix *m); void gretl_matrix_xtr_symmetric (gretl_matrix *m); void gretl_matrix_set_equals_tolerance (double tol); void gretl_matrix_unset_equals_tolerance (void); gretl_matrix * gretl_matrix_alloc (int rows, int cols); gretl_matrix * gretl_matrix_reuse (gretl_matrix *m, int rows, int cols); int gretl_matrix_realloc (gretl_matrix *m, int rows, int cols); void gretl_matrix_block_destroy (gretl_matrix_block *B); gretl_matrix_block * gretl_matrix_block_new (gretl_matrix **pm, ...); gretl_matrix * gretl_identity_matrix_new (int n); gretl_matrix * gretl_DW_matrix_new (int n); gretl_matrix * gretl_zero_matrix_new (int r, int c); gretl_matrix * gretl_unit_matrix_new (int r, int c); gretl_matrix * gretl_null_matrix_new (void); gretl_matrix * gretl_matrix_seq (int start, int end, int step, int *err); gretl_matrix * gretl_matrix_copy (const gretl_matrix *m); int gretl_matrix_copy_row (gretl_matrix *dest, int di, const gretl_matrix *src, int si); int gretl_matrix_inscribe_I (gretl_matrix *m, int row, int col, int n); gretl_matrix * gretl_matrix_copy_transpose (const gretl_matrix *m); gretl_matrix * gretl_matrix_reverse_rows (const gretl_matrix *m); gretl_matrix * gretl_matrix_reverse_cols (const gretl_matrix *m); gretl_matrix * gretl_matrix_get_diagonal (const gretl_matrix *m, int *err); double gretl_matrix_trace (const gretl_matrix *m, int *err); int gretl_matrix_random_fill (gretl_matrix *m, int dist); gretl_matrix * gretl_random_matrix_new (int r, int c, int dist); gretl_matrix * gretl_matrix_resample (const gretl_matrix *m, int *err); gretl_matrix * gretl_matrix_block_resample (const gretl_matrix *m, int blocklen, int *err); double gretl_vector_mean (const gretl_vector *v); double gretl_vector_variance (const gretl_vector *v); void gretl_matrix_zero (gretl_matrix *m); int gretl_matrix_zero_upper (gretl_matrix *m); int gretl_matrix_zero_lower (gretl_matrix *m); void gretl_matrix_fill (gretl_matrix *m, double x); void gretl_matrix_multiply_by_scalar (gretl_matrix *m, double x); int gretl_matrix_divide_by_scalar (gretl_matrix *m, double x); void gretl_matrix_switch_sign (gretl_matrix *m); gretl_matrix * gretl_matrix_dot_op (const gretl_matrix *a, const gretl_matrix *b, int op, int *err); gretl_matrix * gretl_matrix_complex_multiply (const gretl_matrix *a, const gretl_matrix *b, int *err); gretl_matrix * gretl_matrix_complex_divide (const gretl_matrix *a, const gretl_matrix *b, int *err); gretl_matrix * gretl_matrix_exp (const gretl_matrix *m, int *err); gretl_matrix * gretl_matrix_polroots (const gretl_matrix *a, int *err); void gretl_matrix_raise (gretl_matrix *m, double x); void gretl_matrix_free (gretl_matrix *m); double * gretl_matrix_steal_data (gretl_matrix *m); int gretl_vector_copy_values (gretl_vector *targ, const gretl_vector *src); int gretl_matrix_copy_values (gretl_matrix *targ, const gretl_matrix *src); int gretl_matrix_copy_values_shaped (gretl_matrix *targ, const gretl_matrix *src); int gretl_matrix_add_to (gretl_matrix *targ, const gretl_matrix *src); int gretl_matrix_add_transpose_to (gretl_matrix *targ, const gretl_matrix *src); int gretl_matrix_subtract_from (gretl_matrix *targ, const gretl_matrix *src); int gretl_matrix_subtract_reversed (const gretl_matrix *a, gretl_matrix *b); int gretl_matrix_I_minus (gretl_matrix *m); int gretl_matrix_transpose_in_place (gretl_matrix *m); int gretl_matrix_transpose (gretl_matrix *targ, const gretl_matrix *src); int gretl_square_matrix_transpose (gretl_matrix *m); int gretl_matrix_add_self_transpose (gretl_matrix *m); int gretl_matrix_vectorize (gretl_matrix *targ, const gretl_matrix *src); int gretl_matrix_unvectorize (gretl_matrix *targ, const gretl_matrix *src); int gretl_matrix_vectorize_h (gretl_matrix *targ, const gretl_matrix *src); int gretl_matrix_unvectorize_h (gretl_matrix *targ, const gretl_matrix *src); int gretl_matrix_inscribe_matrix (gretl_matrix *targ, const gretl_matrix *src, int row, int col, GretlMatrixMod mod); int gretl_matrix_extract_matrix (gretl_matrix *targ, const gretl_matrix *src, int row, int col, GretlMatrixMod mod); int gretl_matrix_multiply_mod (const gretl_matrix *a, GretlMatrixMod amod, const gretl_matrix *b, GretlMatrixMod bmod, gretl_matrix *c, GretlMatrixMod cmod); int gretl_matrix_multiply (const gretl_matrix *a, const gretl_matrix *b, gretl_matrix *c); gretl_matrix * gretl_matrix_multiply_new (const gretl_matrix *a, const gretl_matrix *b, int *err); int gretl_matrix_kronecker_product (const gretl_matrix *A, const gretl_matrix *B, gretl_matrix *K); gretl_matrix * gretl_matrix_kronecker_product_new (const gretl_matrix *A, const gretl_matrix *B, int *err); int gretl_matrix_I_kronecker (int p, const gretl_matrix *B, gretl_matrix *K); gretl_matrix * gretl_matrix_I_kronecker_new (int p, const gretl_matrix *B, int *err); int gretl_matrix_kronecker_I (const gretl_matrix *A, int r, gretl_matrix *K); gretl_matrix * gretl_matrix_kronecker_I_new (const gretl_matrix *A, int r, int *err); gretl_matrix * gretl_matrix_pow (const gretl_matrix *A, int s, int *err); double gretl_matrix_dot_product (const gretl_matrix *a, GretlMatrixMod amod, const gretl_matrix *b, GretlMatrixMod bmod, int *errp); double gretl_vector_dot_product (const gretl_vector *a, const gretl_vector *b, int *errp); gretl_matrix * gretl_matrix_row_sum (const gretl_matrix *m, int *err); gretl_matrix * gretl_matrix_column_sum (const gretl_matrix *m, int *err); gretl_matrix * gretl_matrix_row_mean (const gretl_matrix *m, int *err); gretl_matrix * gretl_matrix_column_mean (const gretl_matrix *m, int *err); gretl_matrix * gretl_matrix_column_sd (const gretl_matrix *m, int *err); double gretl_matrix_row_i_mean (const gretl_matrix *m, int row); double gretl_matrix_column_j_mean (const gretl_matrix *m, int col); void gretl_matrix_demean_by_row (gretl_matrix *m); void gretl_matrix_demean_by_column (gretl_matrix *m); gretl_matrix * gretl_matrix_vcv (gretl_matrix *m); gretl_matrix * gretl_matrix_quantiles (const gretl_matrix *m, double p, int *err); double gretl_matrix_determinant (gretl_matrix *a, int *err); double gretl_matrix_log_determinant (gretl_matrix *a, int *err); double gretl_matrix_log_abs_determinant (gretl_matrix *a, int *err); double gretl_vcv_log_determinant (const gretl_matrix *m); double gretl_matrix_one_norm (const gretl_matrix *m); double gretl_matrix_infinity_norm (const gretl_matrix *m); int gretl_LU_solve (gretl_matrix *a, gretl_matrix *b); int gretl_cholesky_decomp_solve (gretl_matrix *a, gretl_matrix *b); int gretl_cholesky_solve (const gretl_matrix *a, gretl_vector *b); gretl_vector * gretl_toeplitz_solve (const gretl_vector *c, const gretl_vector *r, const gretl_vector *b, int *err); gretl_matrix * gretl_matrix_XTX_new (const gretl_matrix *X); int gretl_inverse_from_cholesky_decomp (gretl_matrix *targ, const gretl_matrix *src); int gretl_invert_general_matrix (gretl_matrix *a); int gretl_invert_symmetric_indef_matrix (gretl_matrix *a); int gretl_invert_symmetric_matrix (gretl_matrix *a); int gretl_invert_symmetric_matrix2 (gretl_matrix *a, double *ldet); int gretl_invert_packed_symmetric_matrix (gretl_matrix *v); int gretl_invert_triangular_matrix (gretl_matrix *a, char uplo); int gretl_invert_diagonal_matrix (gretl_matrix *a); int gretl_invert_matrix (gretl_matrix *a); int gretl_matrix_moore_penrose (gretl_matrix *A); int gretl_SVD_invert_matrix (gretl_matrix *a); int gretl_invpd (gretl_matrix *a); int gretl_maybe_invpd (gretl_matrix *a); int gretl_matrix_SVD (const gretl_matrix *a, gretl_matrix **pu, gretl_vector **ps, gretl_matrix **pvt); double gretl_symmetric_matrix_rcond (const gretl_matrix *m, int *err); double gretl_matrix_rcond (const gretl_matrix *m, int *err); int gretl_general_eigen_sort (gretl_matrix *evals, gretl_matrix *evecs, int rank); int gretl_symmetric_eigen_sort (gretl_matrix *evals, gretl_matrix *evecs, int rank); gretl_matrix * gretl_general_matrix_eigenvals (gretl_matrix *m, int eigenvecs, int *err); gretl_matrix * gretl_symmetric_matrix_eigenvals (gretl_matrix *m, int eigenvecs, int *err); gretl_matrix * gretl_gensymm_eigenvals (const gretl_matrix *A, const gretl_matrix *B, gretl_matrix *V, int *err); gretl_matrix * gretl_matrix_right_nullspace (const gretl_matrix *M, int *err); gretl_matrix * gretl_matrix_left_nullspace (const gretl_matrix *M, GretlMatrixMod mod, int *err); gretl_matrix * gretl_matrix_row_concat (const gretl_matrix *a, const gretl_matrix *b, int *err); gretl_matrix * gretl_matrix_col_concat (const gretl_matrix *a, const gretl_matrix *b, int *err); int gretl_matrix_inplace_colcat (gretl_matrix *a, const gretl_matrix *b, const char *mask); gretl_matrix * gretl_matrix_cumcol (const gretl_matrix *m, int *err); gretl_matrix * gretl_matrix_diffcol (const gretl_matrix *m, double missval, int *err); gretl_matrix * gretl_matrix_lag (const gretl_matrix *m, const gretl_vector *k, double missval); int gretl_matrix_inplace_lag (gretl_matrix *targ, const gretl_matrix *src, int k); int gretl_matrix_cholesky_decomp (gretl_matrix *a); int gretl_matrix_psd_root (gretl_matrix *a); int gretl_matrix_QR_decomp (gretl_matrix *M, gretl_matrix *R); int gretl_check_QR_rank (const gretl_matrix *R, int *err, double *rcnd); int gretl_matrix_rank (const gretl_matrix *a, int *err); int gretl_matrix_ols (const gretl_vector *y, const gretl_matrix *X, gretl_vector *b, gretl_matrix *vcv, gretl_vector *uhat, double *s2); int gretl_matrix_multi_ols (const gretl_matrix *Y, const gretl_matrix *X, gretl_matrix *B, gretl_matrix *E, gretl_matrix **XTXi); int gretl_matrix_multi_SVD_ols (const gretl_matrix *Y, const gretl_matrix *X, gretl_matrix *B, gretl_matrix *E, gretl_matrix **XTXi); 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); double gretl_matrix_r_squared (const gretl_matrix *y, const gretl_matrix *X, const gretl_matrix *b, int *err); 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); 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); int gretl_matrix_SVD_ols (const gretl_vector *y, const gretl_matrix *X, gretl_vector *b, gretl_matrix *vcv, gretl_vector *uhat, double *s2); int gretl_matrix_qform (const gretl_matrix *A, GretlMatrixMod amod, const gretl_matrix *X, gretl_matrix *C, GretlMatrixMod cmod); double gretl_scalar_qform (const gretl_vector *b, const gretl_matrix *X, int *errp); int gretl_matrix_columnwise_product (const gretl_matrix *A, const gretl_matrix *B, gretl_matrix *C); int gretl_matrix_diagonal_sandwich (const gretl_vector *d, const gretl_matrix *X, gretl_matrix *DXD); void gretl_matrix_set_t1 (gretl_matrix *m, int t); void gretl_matrix_set_t2 (gretl_matrix *m, int t); int gretl_matrix_get_t1 (const gretl_matrix *m); int gretl_matrix_get_t2 (const gretl_matrix *m); int gretl_is_identity_matrix (const gretl_matrix *m); int gretl_is_zero_matrix (const gretl_matrix *m); int gretl_matrix_get_structure (const gretl_matrix *m); int gretl_matrices_are_equal (const gretl_matrix *a, const gretl_matrix *b, int *err); gretl_matrix * gretl_covariance_matrix (const gretl_matrix *m, int corr, int *errp); gretl_matrix ** gretl_matrix_array_new (int n); gretl_matrix ** gretl_matrix_array_new_with_size (int n, int rows, int cols); void gretl_matrix_array_free (gretl_matrix **A, int n); gretl_matrix * gretl_matrix_values (const double *x, int n, int *err); gretl_matrix * gretl_matrix_shape (const gretl_matrix *A, int r, int c); gretl_matrix * gretl_matrix_trim_rows (const gretl_matrix *A, int ttop, int tbot, int *err); gretl_matrix * gretl_matrix_minmax (const gretl_matrix *A, int mm, int rc, int idx, int *err); gretl_matrix * gretl_matrix_pca (const gretl_matrix *X, int p, int *err); gretl_matrix * gretl_matrix_xtab (int t1, int t2, const double *x, const double *y, int *err); gretl_matrix * matrix_matrix_xtab (const gretl_matrix *x, const gretl_matrix *y, int *err); gretl_matrix * gretl_matrix_bool_sel (const gretl_matrix *A, const gretl_matrix *sel, int rowsel, int *err); gretl_matrix * gretl_matrix_sort_by_column (const gretl_matrix *m, int k, int *err); gretl_matrix * gretl_matrix_covariogram (const gretl_matrix *X, const gretl_matrix *u, const gretl_matrix *w, int p, int *err); void gretl_matrix_transcribe_obs_info (gretl_matrix *targ, const gretl_matrix *src); void lapack_mem_free (void); void set_blas_nmk_min (int n); int get_blas_nmk_min (void);
typedef enum {
GRETL_MOD_NONE = 0,
GRETL_MOD_TRANSPOSE,
GRETL_MOD_SQUARE,
GRETL_MOD_CUMULATE,
GRETL_MOD_DECREMENT
} GretlMatrixMod;
typedef enum {
GRETL_MATRIX_SQUARE = 1,
GRETL_MATRIX_LOWER_TRIANGULAR,
GRETL_MATRIX_UPPER_TRIANGULAR,
GRETL_MATRIX_SYMMETRIC,
GRETL_MATRIX_DIAGONAL,
GRETL_MATRIX_IDENTITY,
GRETL_MATRIX_SCALAR,
} GretlMatrixStructure;
#define gretl_matrix_get(m,i,j) (m->val[(j)*m->rows+(i)])
m : |
matrix. |
i : |
row. |
j : |
column. |
| Returns : | the i, j element of m.
|
#define gretl_vector_get(v,i) (v->val[i])
v : |
vector. |
i : |
index. |
| Returns : | element i of v.
|
#define gretl_matrix_set(m,i,j,x) ((m)->val[(j)*(m)->rows+(i)]=x)
Sets the i, j element of m to x.
m : |
matrix. |
i : |
row. |
j : |
column. |
x : |
value to set. |
| Returns : |
#define gretl_vector_set(v,i,x) ((v)->val[i]=x)
Sets element i of v to x.
v : |
vector. |
i : |
index. |
x : |
value to set. |
| Returns : |
#define gretl_column_vector_alloc(i) gretl_matrix_alloc((i),1)
i : |
number of rows. |
#define gretl_vector_free(v) gretl_matrix_free(v)
Frees the vector v and its associated storage.
v : |
gretl_vector to free.
|
#define gretl_is_null_matrix(m) (m == NULL || m->rows == 0 || m->cols == 0)
m : |
void gretl_matrix_print (const gretl_matrix *m, const char *msg);
Prints the given matrix to stderr.
m : |
matrix. |
msg : |
message to print with matrix, or NULL.
|
int gretl_matrix_is_symmetric (const gretl_matrix *m);
m : |
gretl_matrix. |
| Returns : | 1 if m is symmetric (with a small relative tolerance
for asymmetry), otherwise 0.
|
int gretl_matrix_is_idempotent (const gretl_matrix *m);
m : |
gretl_matrix. |
| Returns : | 1 if m is idempotent, otherwise 0.
|
void gretl_matrix_xtr_symmetric (gretl_matrix *m);
Computes the symmetric part of m by averaging its off-diagonal
elements.
m : |
gretl_matrix. |
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.5e-12.
tol : |
tolerance value. |
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 * gretl_matrix_alloc (int rows, int cols);
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_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 the when the matrix is to be freed, gretl_matrix_free()
should be applied only once.
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. |
int gretl_matrix_realloc (gretl_matrix *m, int rows, int cols);
Reallocates the storage in m to the specified dimensions.
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_block * gretl_matrix_block_new (gretl_matrix **pm, ...);
pm : |
|
... : |
|
| Returns : |
gretl_matrix * gretl_identity_matrix_new (int n);
n : |
desired number of rows and columns in the matrix. |
| Returns : | pointer to a newly allocated identity matrix, or NULL
on failure.
|
gretl_matrix * gretl_DW_matrix_new (int n);
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_matrix * gretl_zero_matrix_new (int r, int c);
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_matrix * gretl_unit_matrix_new (int r, int c);
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_matrix * gretl_null_matrix_new (void);
| Returns : | pointer to a newly allocated null matrix
(for use in declaration of a variable as a matrix),
NULL on failure.
|
gretl_matrix * gretl_matrix_seq (int start, int end, int step, int *err);
start : |
first element. |
end : |
last element. |
step : |
positive integer step. |
err : |
location to recieve error code. |
| Returns : | pointer to a row vector, containing the numbers from
start to end, in decreasing order if start > end --
or NULL on failure.
|
gretl_matrix * gretl_matrix_copy (const gretl_matrix *m);
m : |
source matrix to be copied. |
| Returns : | an allocated copy of matrix m, or NULL on failure.
|
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.
dest : |
destination matrix. |
di : |
row to copy into. |
src : |
source matrix. |
si : |
row to copy from. |
| Returns : | 0 on success, non-zero on failure. |
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).
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 * gretl_matrix_copy_transpose (const gretl_matrix *m);
m : |
source matrix to be copied. |
| Returns : | an allocated copy of the tranpose of m, or NULL on failure.
|
gretl_matrix * gretl_matrix_reverse_rows (const gretl_matrix *m);
m : |
source matrix whose rows are to be reversed. |
| Returns : | a matrix with the same rows as m, last to first.
|
gretl_matrix * gretl_matrix_reverse_cols (const gretl_matrix *m);
m : |
source matrix whose columns are to be reversed. |
| Returns : | a matrix with the same columns as m, last to first.
|
gretl_matrix * gretl_matrix_get_diagonal (const gretl_matrix *m, int *err);
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.
|
double gretl_matrix_trace (const gretl_matrix *m, int *err);
m : |
square input matrix. |
err : |
location to receive error code. |
| Returns : | the trace (sum of diagonal elements) of m, if
m is square, otherwise NADBL.
|
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.
m : |
input matrix. |
dist : |
either D_UNIFORM or D_NORMAL.
|
| Returns : | 0 on success, 1 on failure. |
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.
r : |
number of rows. |
c : |
number of columns. |
dist : |
either D_UNIFORM or D_NORMAL.
|
| Returns : | allocated matrix or NULL on failure.
|
gretl_matrix * gretl_matrix_resample (const gretl_matrix *m, int *err);
m : |
input matrix. |
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 * gretl_matrix_block_resample (const gretl_matrix *m, int blocklen, int *err);
m : |
input matrix. |
blocklen : |
length of moving blocks. |
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.
|
double gretl_vector_mean (const gretl_vector *v);
v : |
input vector. |
| Returns : | the arithmetic mean of the elements of v, or
NADBL on failure.
|
double gretl_vector_variance (const gretl_vector *v);
v : |
input vector. |
| Returns : | the variance of the elements of v, or
NADBL on failure.
|
void gretl_matrix_zero (gretl_matrix *m);
Sets all elements of m to zero.
m : |
matrix to be set to zero. |
int gretl_matrix_zero_upper (gretl_matrix *m);
Sets the elements of m outside of the lower triangle to zero.
m : |
square matrix to operate on. |
| Returns : | 0 on success, non-zero error code otherwise. |
int gretl_matrix_zero_lower (gretl_matrix *m);
Sets the elements of m outside of the upper triangle to zero.
m : |
square matrix to operate on. |
| Returns : | 0 on success, non-zero error code otherwise. |
void gretl_matrix_fill (gretl_matrix *m, double x);
Sets all entries in m to the value x.
m : |
matrix to fill. |
x : |
value with which to fill. |
void gretl_matrix_multiply_by_scalar (gretl_matrix *m, double x);
Multiplies all elements of m by x.
m : |
matrix to operate on. |
x : |
scalar by which to multiply. |
int gretl_matrix_divide_by_scalar (gretl_matrix *m, double x);
Divides all elements of m by x.
m : |
matrix to operate on. |
x : |
scalar by which to divide. |
| Returns : | 0 on success, 1 if x = 0. |
void gretl_matrix_switch_sign (gretl_matrix *m);
Changes the sign of each element of m.
m : |
matrix to operate on. |
gretl_matrix * gretl_matrix_dot_op (const gretl_matrix *a, const gretl_matrix *b, int op, int *err);
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).
|
gretl_matrix * gretl_matrix_complex_multiply (const gretl_matrix *a, const gretl_matrix *b, 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.
a : |
m x (1 or 2) matrix. |
b : |
m x (1 or 2) matrix. |
err : |
location to receive error code. |
| Returns : | an m x 2 matrix with the result of the multiplication
of the two vectors of complex numbers. If both a and b have no
imaginary part, the return value will be m x 1. Or NULL on
failure.
|
gretl_matrix * gretl_matrix_complex_divide (const gretl_matrix *a, const gretl_matrix *b, 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.
a : |
m x (1 or 2) matrix. |
b : |
m x (1 or 2) matrix. |
err : |
location to receive error code. |
| Returns : | an m x 2 matrix with the result of the division
of the two vectors of complex numbers. If both a and b have no
imaginary part, the return value will be m x 1. Or NULL on
failure.
|
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.
m : |
square matrix to operate on. |
err : |
location to receive error code. |
| Returns : | the exponential, or NULL on failure.
|
gretl_matrix * gretl_matrix_polroots (const gretl_matrix *a, 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.
a : |
vector of coefficients. |
err : |
location to receive error code. |
| Returns : | a p-vector if all the roots are real, otherwise a
p x 2 matrix with the real parts in the first column and
the imaginary parts in the second. Or NULL on failure.
|
void gretl_matrix_raise (gretl_matrix *m, double x);
Raises each element of m to the power x.
m : |
matrix to operate on. |
x : |
exponent. |
void gretl_matrix_free (gretl_matrix *m);
Frees the allocated storage in m, then m itself.
m : |
matrix to be freed. |
double * gretl_matrix_steal_data (gretl_matrix *m);
"Steals" the allocated data from m, which is left with a
NULL data pointer.
m : |
matrix to operate on. |
| Returns : | a pointer to the "stolen" data. |
int gretl_vector_copy_values (gretl_vector *targ, const gretl_vector *src);
Copies the elements of src into the corresponding elements
of targ.
targ : |
target vector. |
src : |
source vector. |
| Returns : | 0 on successful completion, or E_NONCONF if the
two vectors are not of the same length.
|
int gretl_matrix_copy_values (gretl_matrix *targ, const gretl_matrix *src);
Copies the elements of src into the corresponding elements
of targ.
targ : |
target matrix. |
src : |
source matrix. |
| Returns : | 0 on successful completion, or
E_NONCONF if the two matrices are not
conformable for the operation.
|
int gretl_matrix_copy_values_shaped (gretl_matrix *targ, const gretl_matrix *src);
Copies the elements of src into targ, column
by column.
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.
|
int gretl_matrix_add_to (gretl_matrix *targ, const gretl_matrix *src);
Adds the elements of src to the corresponding elements
of targ.
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.
|
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.
targ : |
target matrix. |
src : |
source matrix. |
| Returns : | 0 on successful completion, or E_NONCONF if the
two matrices are not conformable for the operation.
|
int gretl_matrix_subtract_from (gretl_matrix *targ, const gretl_matrix *src);
Subtracts the elements of src from the corresponding elements
of targ.
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.
|
int gretl_matrix_subtract_reversed (const gretl_matrix *a, gretl_matrix *b);
Operates on b such that b_{ij} = a_{ij} - b_{ij}.
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.
|
int gretl_matrix_I_minus (gretl_matrix *m);
Rewrites m as (I - m), where I is the n x n identity
matrix.
m : |
original square matrix, n x n. |
| Returns : | 0 on successful completion, or E_NONCONF if m is
not square.
|
int gretl_matrix_transpose_in_place (gretl_matrix *m);
Tranposes m in place.
m : |
matrix to transpose. |
| Returns : | 0 on success, non-zero error code otherwise. |
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.
targ : |
target matrix. |
src : |
source matrix. |
| Returns : | 0 on success, non-zero error code otherwise. |
int gretl_square_matrix_transpose (gretl_matrix *m);
Transposes the matrix m.
m : |
square matrix to operate on. |
| Returns : | 0 on success, non-zero error code otherwise. |
int gretl_matrix_add_self_transpose (gretl_matrix *m);
Adds the transpose of m to m itself, yielding a symmetric
matrix.
m : |
(square) matrix to operate on. |
| Returns : | 0 on successful completion, or 1 if the source matrix is not square. |
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.
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.
|
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).
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.
|
int gretl_matrix_vectorize_h (gretl_matrix *targ, const gretl_matrix *src);
Writes into targ vech(src), that is, a column 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.
targ : |
target vector, (m * (m+1)/2) x 1. |
src : |
source square matrix, m x m. |
| Returns : | 0 on successful completion, or E_NONCONF if
targ is not correctly dimensioned.
|
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.
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.
|
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.
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.
|
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.
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.
|
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.
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.
|
| Returns : | 0 on success; non-zero error code on failure. |
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.
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 * 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.
a : |
left-hand matrix. |
b : |
right-hand matrix. |
err : |
location for error code. |
| Returns : | matrix product on success, or NULL on failure.
|
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.
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 * gretl_matrix_kronecker_product_new (const gretl_matrix *A, const gretl_matrix *B, int *err);
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.
|
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.
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 * 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.
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.
|
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.
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 * 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.
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 * gretl_matrix_pow (const gretl_matrix *A, int s, int *err);
Calculates the matrix A^k using Golub and Van Loan's Algorithm 11.2.2 ("Binary Powering").
A : |
square source matrix. |
s : |
exponent >= 0. |
err : |
location to receive error code. |
| Returns : | allocated matrix, or NULL on failure.
|
double gretl_matrix_dot_product (const gretl_matrix *a, GretlMatrixMod amod, const gretl_matrix *b, GretlMatrixMod bmod, int *errp);
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.
|
errp : |
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.
|
double gretl_vector_dot_product (const gretl_vector *a, const gretl_vector *b, int *errp);
gretl_matrix * gretl_matrix_row_sum (const gretl_matrix *m, int *err);
m : |
source matrix. |
err : |
|
| Returns : | a column vector containing the sums of
the rows of m, or NULL on failure.
|
gretl_matrix * gretl_matrix_column_sum (const gretl_matrix *m, int *err);
m : |
source matrix. |
err : |
|
| Returns : | a row vector containing the sums of
the columns of m, or NULL on failure.
|
gretl_matrix * gretl_matrix_row_mean (const gretl_matrix *m, int *err);
m : |
source matrix. |
err : |
location to receive error code. |
| Returns : | a column vector containing the means of
the rows of m, or NULL on failure.
|
gretl_matrix * gretl_matrix_column_mean (const gretl_matrix *m, int *err);
m : |
source matrix. |
err : |
location to receive error code. |
| Returns : | a row vector containing the means of
the columns of m, or NULL on failure.
|
gretl_matrix * gretl_matrix_column_sd (const gretl_matrix *m, int *err);
m : |
source matrix. |
err : |
location to receive error code. |
| Returns : | a row vector containing the standard deviations of
the columns of m (without a degrees of freedom correction),
or NULL on failure.
|
double gretl_matrix_row_i_mean (const gretl_matrix *m, int row);
m : |
source matrix. |
row : |
zero-based index of row. |
| Returns : | the mean of the elements in row row of m,
or NADBL if the row is out of bounds.
|
double gretl_matrix_column_j_mean (const gretl_matrix *m, int col);
m : |
source matrix. |
col : |
zero-based index of column. |
| Returns : | the mean of the elements in column col of m,
or NADBL if the column is out of bounds.
|
void gretl_matrix_demean_by_row (gretl_matrix *m);
For each row of m, subtracts the row mean from each
element on the row.
m : |
matrix on which to operate. |
void gretl_matrix_demean_by_column (gretl_matrix *m);
For each column of m, subtracts the column mean from each
element in the column.
m : |
matrix on which to operate. |
gretl_matrix * gretl_matrix_vcv (gretl_matrix *m);
Forms a variance-covariance matrix based on m, thus:
(1) subtract the column means from the column elements of m;
(2) multiply m-transpose into m; and
(3) divide the elements of the product by the number of rows
in m.
m : |
source matrix (expected to have rows >= cols). |
| Returns : | the allocated variance-covariance matrix, or NULL
on failure. Note that on return the column means have
been subtracted from m.
|
gretl_matrix * gretl_matrix_quantiles (const gretl_matrix *m, double p, int *err);
m : |
matrix on which to operate. |
p : |
probability. |
err : |
location to receive error code. |
| Returns : | a row vector containing the p quantiles
of the columns of m, or NULL on failure.
|
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.
a : |
gretl_matrix. |
err : |
location to receive error code. |
| Returns : | the determinant, or NABDL on failure. |
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.
a : |
gretl_matrix. |
err : |
location to receive error code. |
| Returns : | the determinant, or NABDL on failure. |
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.
a : |
gretl_matrix. |
err : |
location to receive error code. |
| Returns : | the determinant, or NABDL on failure. |
double gretl_vcv_log_determinant (const gretl_matrix *m);
Compute the log determinant of the symmetric positive-definite
matrix m using Cholesky decomposition.
m : |
gretl_matrix. |
| Returns : | the log determinant, or NABDL on failure. |
double gretl_matrix_one_norm (const gretl_matrix *m);
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).
|
double gretl_matrix_infinity_norm (const gretl_matrix *m);
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).
|
int gretl_LU_solve (gretl_matrix *a, gretl_matrix *b);
Solves ax = b for the unknown x, using LU decomposition.
On exit, b is replaced by the solution and a is replaced
by its LU decomposition.
a : |
gretl_matrix. |
b : |
gretl_matrix. |
| Returns : | 0 on successful completion, non-zero code on error. |
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.
a : |
symmetric positive-definite matrix. |
b : |
vector 'x'. |
| Returns : | 0 on successful completion, or non-zero code on error. |
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.
a : |
Cholesky-decomposed symmetric positive-definite matrix. |
b : |
vector 'x'. |
| Returns : | 0 on successful completion, or non-zero code on error. |
gretl_vector * gretl_toeplitz_solve (const gretl_vector *c, const gretl_vector *r, const gretl_vector *b, 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].
c : |
Toeplitz column. |
r : |
Toeplitz row. |
b : |
right-hand side vector. |
err : |
error code. |
| Returns : | a newly allocated vector, containing the solution, x. |
gretl_matrix * gretl_matrix_XTX_new (const gretl_matrix *X);
X : |
matrix to process. |
| Returns : | a newly allocated matrix containing X'X, or
NULL on error.
|
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.
targ : |
matrix to hold inverse. |
src : |
Cholesky-decomposed matrix. |
| Returns : | 0 on success; non-zero error code on failure. |
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.
a : |
matrix to invert. |
| Returns : | 0 on success; non-zero error code on failure. |
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.
a : |
matrix to invert. |
| Returns : | 0 on success; non-zero error code on failure. |
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.
a : |
matrix to invert. |
| Returns : | 0 on success; non-zero error code on failure. |
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.
a : |
matrix to invert. |
ldet : |
location to recieve log determinant, or NULL.
|
| Returns : | 0 on success; non-zero error code on failure. |
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.
v : |
symmetric matrix in vech form (lower triangle packed as a column vector). |
| Returns : | 0 on success; non-zero error code on failure. |
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.
a : |
triangular matrix to invert. |
uplo : |
'L' for lower triangular a, 'U' for upper.
|
| Returns : | 0 on success; non-zero error code on failure. |
int gretl_invert_diagonal_matrix (gretl_matrix *a);
Computes the inverse of a diagonal matrix.
On exit a is overwritten with the inverse.
a : |
matrix to invert. |
| Returns : | 0 on success; non-zero error code on failure. |
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.
a : |
matrix to invert. |
| Returns : | 0 on success; non-zero error code on failure. |
int gretl_matrix_moore_penrose (gretl_matrix *A);
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.
A : |
m x n matrix. |
| Returns : | 0 on success; non-zero error code on failure. |
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.
a : |
n x n matrix to invert. |
| Returns : | 0 on success; non-zero error code on failure. |
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.
a : |
matrix to invert. |
| Returns : | 0 on success; non-zero error code on failure. |
int gretl_maybe_invpd (gretl_matrix *a);
Attempts to computes the inverse of a matrix which may be
positive definite. 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. Unlike gretl_invpd() this function
does not sump error messages to stderr in case the matrix
is not in fact positive definite.
a : |
matrix to invert. |
| Returns : | 0 on success; non-zero error code on failure. |
int gretl_matrix_SVD (const gretl_matrix *a, gretl_matrix **pu, gretl_vector **ps, gretl_matrix **pvt);
Computes SVD factorization of a general matrix using the lapack
function dgesvd. A = u * diag(s) * vt.
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.
m : |
matrix to examine. |
err : |
location to receive error code. |
| Returns : | the estimate, or NADBL on failure to allocate memory. |
double gretl_matrix_rcond (const gretl_matrix *m, int *err);
Estimates the reciprocal condition number of the real
matrix m (in the 1-norm).
m : |
matrix to examine. |
err : |
location to receive error code. |
| Returns : | the estimate, or NADBL on failure to allocate memory. |
int gretl_general_eigen_sort (gretl_matrix *evals, gretl_matrix *evecs, int rank);
Sorts the real components of 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.
evals : |
array of eigenvalues, from general (not necessarily symmetric) matrix. |
evecs : |
matrix of eigenvectors. |
rank : |
desired number of columns in output. |
| Returns : | 0 on success; non-zero error code on failure. |
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.
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_matrix * gretl_general_matrix_eigenvals (gretl_matrix *m, int eigenvecs, int *err);
Computes the eigenvalues of the general matrix m.
If eigenvecs is non-zero, also compute the right
eigenvectors of m, which are stored in m. Uses the lapack
function dgeev.
m : |
square matrix on which to operate. |
eigenvecs : |
non-zero to calculate eigenvectors, 0 to omit. |
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 contains the real parts of
the eigenvalues of m, and the second holds the
imaginary parts.
|
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 dsyev.
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, or NULL
on failure.
|
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.
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).
M : |
matrix to operate on. |
err : |
location to receive error code. |
| Returns : | the allocated matrix R, or NULL on failure.
|
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).
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 * gretl_matrix_row_concat (const gretl_matrix *a, const gretl_matrix *b, int *err);
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 * gretl_matrix_col_concat (const gretl_matrix *a, const gretl_matrix *b, int *err);
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.
|
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.
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 * gretl_matrix_cumcol (const gretl_matrix *m, int *err);
m : |
source matrix. |
err : |
error code. |
| Returns : | a matrix of the same dimensions as m, containing
the cumulated columns of m.
|
gretl_matrix * gretl_matrix_diffcol (const gretl_matrix *m, double missval, int *err);
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 * gretl_matrix_lag (const gretl_matrix *m, const gretl_vector *k, double missval);
m : |
source matrix. |
k : |
vector of lag orders (> 0 for lags, < 0 for leads). |
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.
|
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.
targ : |
target matrix. |
src : |
source matrix. |
k : |
lag order (> 0 for lags, < 0 for leads). |
| Returns : | 0 on success, non-zero code otherwise. |
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.
a : |
matrix to operate on. |
| Returns : | 0 on success; 1 on failure. |
int gretl_matrix_psd_root (gretl_matrix *a);
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.
a : |
matrix to operate on. |
| Returns : | 0 on success; 1 on failure. |
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.
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. |
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.
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. |
int gretl_matrix_rank (const gretl_matrix *a, int *err);
Computes the rank of a via its SV decomposition.
a : |
matrix to examine. |
err : |
location to receive error code on failure. |
| Returns : | the rank of a, or 0 on failure.
|
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.
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. |
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.
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} on output, or NULL if this is not
needed.
|
| Returns : | 0 on success, non-zero error code on failure. |
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.
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.
|
| Returns : | 0 on success, non-zero error code on failure. |
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.
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. |
double gretl_matrix_r_squared (const gretl_matrix *y, const gretl_matrix *X, const gretl_matrix *b, int *err);
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 one the regression
represented by y, X and b, or NADBL on failure.
|
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 Johnsen 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.
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. |
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 LU factorization,
and puts the coefficient estimates in b. Optionally, calculates the
covariance matrix in vcv.
y : |
dependent variable vector. |
X : |
matrix of independent variables. |
R : |
left-hand restriction matrix, as in Rb = q. |
q : |
right-hand restriction vector. |
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 ro 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. |
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.
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. |
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.
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.
|
| Returns : | 0 on success; non-zero error code on failure. |
double gretl_scalar_qform (const gretl_vector *b, const gretl_matrix *X, int *errp);
Computes the scalar product bXb', or b'Xb if b is a column
vector. The content of errp is set to 0 on success,
or a non-zero code on failure.
b : |
k-vector. |
X : |
symmetric k x k matrix. |
errp : |
pointer to receive error code. |
| Returns : | the scalar product, or NADBL on failure. |
int gretl_matrix_columnwise_product (const gretl_matrix *A, const gretl_matrix *B, gretl_matrix *C);
A : |
|
B : |
|
C : |
|
| Returns : |
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.
d : |
k-vector. |
X : |
k * k matrix. |
DXD : |
target k * k matrix. |
| Returns : | 0 on success, non-zero code on error. |
void gretl_matrix_set_t1 (gretl_matrix *m, int t);
Sets an integer value on the t1 member of the gretl_matrix
(used for internal information).
m : |
matrix to operate on. |
t : |
integer value to set. |
void gretl_matrix_set_t2 (gretl_matrix *m, int t);
Sets an integer value on the t2 member of the gretl_matrix
(used for internal information).
m : |
matrix to operate on. |
t : |
integer value to set. |
int gretl_matrix_get_t1 (const gretl_matrix *m);
m : |
matrix to read from. |
| Returns : | the integer that has been set on the t1 member
of the matrix struct, or zero.
|
int gretl_matrix_get_t2 (const gretl_matrix *m);
m : |
matrix to read from. |
| Returns : | the integer that has been set on the t2 member
of the matrix struct, or zero.
|
int gretl_is_identity_matrix (const gretl_matrix *m);
m : |
matrix to examine. |
| Returns : | 1 if m is an identity matrix, 0 otherwise.
|
int gretl_is_zero_matrix (const gretl_matrix *m);
m : |
matrix to examine. |
| Returns : | 1 if m is a zero matrix, 0 otherwise.
|
int gretl_matrices_are_equal (const gretl_matrix *a, const gretl_matrix *b, int *err);
a : |
first matrix in comparison. |
b : |
second matrix in comparison. |
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_matrix * gretl_covariance_matrix (const gretl_matrix *m, int corr, int *errp);
m : |
(x x n) matrix containing n observations on each of k variables. |
corr : |
flag for computing correlations. |
errp : |
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 ** 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.
n : |
number of matrices. |
| Returns : | pointer on sucess, NULL on failure.
|
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.
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.
|
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().
A : |
dyamically allocated array of gretl matrices. |
n : |
number of matrices in array. |
gretl_matrix * gretl_matrix_values (const double *x, int n, int *err);
x : |
array to process. |
n : |
length of array. |
err : |
location to receive error code. |
| Returns : | an allocated matrix containing the distinct
values in array x, or NULL on failure.
|
gretl_matrix * gretl_matrix_shape (const gretl_matrix *A, int r, int c);
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.
A : |
array to process. |
r : |
rows of target matrix. |
c : |
columns of target matrix. |
| Returns : | the generated matrix, or NULL on failure.
|
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.
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 * gretl_matrix_minmax (const gretl_matrix *A, int mm, int rc, int idx, 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.
A : |
m x n matrix to examine. |
mm : |
0 for minimum, 1 for maximum. |
rc : |
0 for row, 1 for column. |
idx : |
0 for values, 1 for indices. |
err : |
location to receive error code. |
| Returns : | the generated matrix, or NULL on failure.
|
gretl_matrix * gretl_matrix_pca (const gretl_matrix *X, int p, 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.
X : |
T x m data matrix. |
p : |
number of principal components to return: 0 < p <= m, or p = -1 to return components for eigenvalues > 1.0. |
err : |
location to receive error code. |
| Returns : | the generated matrix, or NULL on failure.
|
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.
t1 : |
start |
t2 : |
end |
x : |
data vector |
y : |
data vector |
err : |
error code |
| Returns : | the generated matrix, or NULL on failure.
|
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.
x : |
data vector |
y : |
data vector |
err : |
error code |
| Returns : | the generated matrix, or NULL on failure.
|
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.
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 * 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.
m : |
matrix. |
k : |
column by which to sort. |
err : |
location to receive error code. |
| Returns : | the generated matrix, or NULL on failure.
|
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.
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.)
targ : |
target matrix. |
src : |
source matrix. |