Top |
int gretl_sort_by (const double *x
,const double *y
,double *z
,const DATASET *dset
);
int diff_series (const double *x
,double *y
,int f
,const DATASET *dset
);
Calculates the differenced counterpart to the input
series x
. If f
= F_SDIFF, the seasonal difference is
computed; if f
= F_LDIFF, the log difference, and if
f
= F_DIFF, the ordinary first difference.
int interpolate_series (const double *x
,double *y
,const DATASET *dset
);
Performs linear interpolation of missing values in x
,
writing the result into y
. Panel data are handled:
interpolation is strictly in the time-series dimension.
No extrapolation is performed.
int standardize_series (const double *x
,double *y
,int dfc
,const DATASET *dset
);
By default calculates the standardized counterpart to the input
series x
, but if dfc
< 0 the result is just centered.
int orthdev_series (const double *x
,double *y
,const DATASET *dset
);
Calculates in y
the forward orthogonal deviations of the input
series x
. That is, y[t] is the scaled difference between x[t]
and the mean of the subsequent observations on x.
int block_resample_series (const double *x
,double *y
,int blocklen
,const DATASET *dset
);
int fracdiff_series (const double *x
,double *y
,double d
,int diff
,int obs
,const DATASET *dset
);
Calculates the fractionally differenced or lagged
counterpart to the input series x
. The fractional
difference operator is defined as (1-L)^d, while the
fractional lag operator 1-(1-L)^d.
x |
array of original data. |
|
y |
array into which to write the result. |
|
d |
fraction by which to difference. |
|
diff |
boolean variable 1 for fracdiff, 0 for fraclag |
|
obs |
used for autoreg calculation, -1 if whole series should be calculated otherwise just the observation for obs is calculated |
|
dset |
data set information. |
int boxcox_series (const double *x
,double *y
,double d
,const DATASET *dset
);
Calculates in y
the Box-Cox transformation for the
input series x
.
int filter_series (const double *x
,double *y
,const DATASET *dset
,gretl_matrix *A
,gretl_matrix *C
,double y0
,gretl_matrix *x0
);
Filters x according to y_t = C(L)/A(L) x_t. If the intended
AR order is p, A
should be a vector of length p. If the
intended MA order is q, C
should be vector of length (q+1),
the first entry giving the coefficient at lag 0. However, if
C
is NULL this is taken to mean that the lag-0 MA coefficient
is unity (and all others are zero).
gretl_matrix * filter_matrix (gretl_matrix *X
,gretl_vector *A
,gretl_vector *C
,double y0
,gretl_matrix *x0
,int *err
);
Filters the columns of x according to y_t = C(L)/A(L) x_t. If the
intended AR order is p, A
should be a vector of length p. If the
intended MA order is q, C
should be vector of length (q+1), the
first entry giving the coefficient at lag 0. However, if C
is
NULL this is taken to mean that the lag-0 MA coefficient is unity
(and all others are zero).
int exponential_movavg_series (const double *x
,double *y
,const DATASET *dset
,double d
,int n
,double y0
);
x |
array of original data. |
|
y |
array into which to write the result. |
|
dset |
dataset information. |
|
d |
coefficient on lagged |
|
n |
number of |
|
y0 |
optional initial value for EMA (or NADBL). |
int movavg_series (const double *x
,double *y
,const DATASET *dset
,int k
,int center
);
int seasonally_adjust_series (const double *x
,double *y
,const char *vname
,DATASET *dset
,int tramo
,gretl_bundle *b
,PRN *prn
);
int panel_statistic (const double *x
,double *y
,const DATASET *dset
,int k
,const double *mask
);
Given the data in x
, constructs in y
a series containing
a panel-data statistic.
x |
source data. |
|
y |
target into which to write. |
|
dset |
data set information. |
|
k |
code representing the desired statistic: F_PNOBS, F_PMIN, F_PMAX, F_PSUM, F_PMEAN, F_PXSUM, F_PSD, F_PXNOBS or F_PXSUM. |
|
mask |
either NULL or a series with 0s for observations to be excluded from the calculations, non-zero values at other observations. |
gretl_matrix * panel_shrink (const double *x
,int skip
,const DATASET *dset
,int *err
);
By default, constructs a column vector holding the first
non-missing observation of x
for each panel unit within
the current sample range (hence skipping any units that
have no valid observations). However, if noskip
is non-zero,
the vector contains an NA for units with all missing
values.
int panel_expand (const gretl_matrix *x
,double *y
,gretlopt opt
,const DATASET *dset
);
Constructs a series by repeating the values in x
,
by default across time (in which case the source
vector must have as many values as there are
individuals in the current sample range).
int hp_filter (const double *x
,double *hp
,const DATASET *dset
,double lambda
,gretlopt opt
);
Calculates the "cycle" component of the time series in
array x
, using the Hodrick-Prescott filter. Adapted from the
original FORTRAN code by E. Prescott.
int oshp_filter (const double *x
,double *hp
,const DATASET *dset
,double lambda
,gretlopt opt
);
Calculates the "cycle" component of the time series in
array x
, using the a one-sided Hodrick-Prescott filter.
The implementation uses the Kalman filter.
int bkbp_filter (const double *x
,double *bk
,const DATASET *dset
,int bkl
,int bku
,int k
);
Calculates the Baxter and King bandpass filter. If bku exceeds the number of available observations, then the "low-pass" version will be computed (weights sum to 1). Otherwise, the ordinary bandpass filter will be applied (weights sum to 0).
int butterworth_filter (const double *x
,double *bw
,const DATASET *dset
,int n
,double cutoff
);
Calculates the Butterworth filter. The code that does this is based on D.S.G. Pollock's IDEOLOG -- see http://www.le.ac.uk/users/dsgp1/
int poly_trend (const double *x
,double *fx
,const DATASET *dset
,int order
);
Calculates a trend via the method of orthogonal polynomials. Based on C code for D.S.G. Pollock's DETREND program.
int weighted_poly_trend (const double *x
,double *fx
,const DATASET *dset
,int order
,gretlopt opt
,double wratio
,double midfrac
);
Calculates a trend via the method of orthogonal polynomials, using the specified weighting scheme. Based on C code for D.S.G. Pollock's DETREND program.
x |
array of original data. |
|
fx |
array into which to write the fitted series. |
|
dset |
data set information. |
|
order |
desired polynomial order. |
|
opt |
weighting option (OPT_Q = quadratic, OPT_B = cosine bell, OPT_C = crenelated). |
|
wratio |
ratio of maximum to minimum weight. |
|
midfrac |
proportion of the data to be treated specially, in the middle. |
void poly_weights (double *w
,int T
,double wmax
,double midfrac
,gretlopt opt
);
Calculates a set of weights; intended for use with polynomial trend fitting.
gretl_matrix * butterworth_gain (int n
,double cutoff
,int hipass
);
int gen_seasonal_dummies (DATASET *dset
,int ref
,int center
);
Adds to the data set (if these variables are not already present) a set of seasonal dummy variables.
int * seasonals_list (DATASET *dset
,int ref
,int center
,int *err
);
Adds to the data set (if these variables are not already present) a set of seasonal dummy variables.
int gen_panel_dummies (DATASET *dset
,gretlopt opt
,PRN *prn
);
Adds to the data set a set of dummy variables corresponding to either the cross-sectional units in a panel, or the time periods.
dset |
dataset struct. |
|
opt |
|
|
prn |
printer for warning, or NULL. |
int gen_unit (DATASET *dset
,int *vnum
);
(For panel data only) adds to the data set an index variable that uniquely identifies the cross-sectional units.
int gen_time (DATASET *dset
,int tm
,int *vnum
);
Generates (and adds to the dataset, if it's not already
present) a time-trend or index variable. This function
is panel-data aware: if the dataset is a panel and
tm
is non-zero, the trend will not simply run
consecutively over the entire range of the data, but
will correctly represent the location in time of
each observation. The index is 1-based.
const double * gretl_plotx (const DATASET *dset
,gretlopt opt
);
Finds or creates a special dummy variable for use on the
x-axis in plotting; this will have the full length of the
data series as given in dset
, and will be appropriately
configured for the data frequency. Do not try to free this
variable.
double * get_fit_or_resid (const MODEL *pmod
,DATASET *dset
,ModelDataIndex idx
,char *vname
,gchar **pdesc
,int *err
);
Creates a full-length array holding the specified model
data, and writes name and description into the vname
and
vlabel
.
int list_linear_combo (double *y
,const int *list
,const gretl_vector *b
,const DATASET *dset
);
gretl_matrix * midas_gradient (int p
,const gretl_matrix *m
,int method
,int *err
);
gretl_matrix * midas_multipliers (gretl_bundle *mb
,int cumulate
,int idx
,int *err
);
int midas_linear_combo (double *y
,const int *list
,const gretl_matrix *theta
,int method
,const DATASET *dset
);
int * vector_to_midas_list (const gretl_matrix *v
,int f_ratio
,const char *prefix
,DATASET *dset
,int *err
);
gretl_matrix * multi_acf (const gretl_matrix *m
,const int *list
,const DATASET *dset
,int p
,int *err
);
gretl_matrix * multi_xcf (const void *px
,int xtype
,const void *py
,int ytype
,const DATASET *dset
,int p
,int *err
);
gretl_matrix * forecast_stats (const double *y
,const double *f
,int t1
,int t2
,int *n_used
,gretlopt opt
,int *err
);
gretl_matrix * matrix_fc_stats (const double *y
,const gretl_matrix *F
,gretlopt opt
,int *err
);
gretl_matrix * duration_func (const double *y
,const double *cens
,int t1
,int t2
,gretlopt opt
,int *err
);
int nadaraya_watson (const double *y
,const double *x
,double h
,DATASET *dset
,int LOO
,double trim
,double *m
);
Implements the Nadaraya-Watson nonparametric estimator for the
conditional mean of y
given x
via the formula
\widehat{m}_h(x)=\frac{\sum_{i=1}^n K_h(x-X_i) Y_i}{\sum_{i=1}^nK_h(x-X_i)}
and computes it for all elements of x
. Note that, in principle,
the kernel K(X_i-X_j) must be computed for every combination of i
and j, but since the function K()
is assumed to be symmetric, we
compute it once to save time.
The scalar h
holds the kernel bandwidth; if LOO
is non-zero the
"leave-one-out" variant of the estimator (essentially a jackknife
estimator; see Pagan and Ullah, Nonparametric Econometrics, p. 119)
is computed.
A rudimentary form of trimming is implemented: the kernel function
is set to 0 when the product of the trim
parameter times the
bandwidth h
exceeds |X_i - X_j|, so as to speed computation and
enhance numerical stability.
int gretl_loess (const double *y
,const double *x
,int poly_order
,double bandwidth
,gretlopt opt
,DATASET *dset
,double *m
);
gretl_matrix * aggregate_by (const double *x
,const double *y
,const int *xlist
,const int *ylist
,const char *fncall
,const DATASET *dset
,int *err
);
Aggregates one or more data series (x) "by" the values of
one or more discrete series (y). In general either x
or
xlist
should be non-NULL, and one of y
or ylist
should
be non-NULL. (If xlist
is non-NULL then x
is ignored,
and similarly for ylist
and y
). For an account of the
matrix that is returned, see the help for gretl's
"aggregate" command.
int fill_day_of_week_array (double *dow
,const double *y
,const double *m
,const double *d
,const DATASET *dset
);
int fill_isoweek_array (double *wknum
,const double *y
,const double *m
,const double *d
,const DATASET *dset
);
gretl_matrix * empirical_cdf (const double *y
,int n
,int *err
);
Calculates the empirical CDF of y
, skipping any missing values.
gretl_matrix * gretl_matrix_vector_stat (const gretl_matrix *m
,GretlVecStat vs
,int rowwise
,int skip_na
,int *err
);
int substitute_values (double *dest
,const double *src
,int n
,const double *v0
,int n0
,const double *v1
,int n1
);
For each of the n
elements in src
, determine if it appears
in v0
: if so, write to targ
the corresponding element of
v1
, otherwise write the src
element to targ
. The method
employed is either "naive" lookup, or if the problem is of
sufficient size a binary tree.
int * list_from_matrix (const gretl_matrix *m
,const char *prefix
,DATASET *dset
,int *err
);
Constructs a series from each column of m
and returns
a list of the series.
gretl_matrix * omega_from_R (const gretl_matrix *R
,int *err
);
Expresses the (n x n) correlation matrix R in spherical coordinates, via m angles, where m = n*(n-1)/2.
gretl_matrix * R_from_omega (const gretl_matrix *omega
,int cholesky
,gretl_matrix *J
,int *err
);
Returns an (n x n) correlation matrix R given its spherical-coordinates representation omega, which should contain m numbers between 0 and M_PI, where m = n*(n-1)/2.
gretl_matrix * felogit_rec_loglik (int t
,int T
,const gretl_matrix *U
,const gretl_matrix *xi
);
t
: int, number of ones
T
: int, number of observations
U
: matrix with the index function
X
: matrix with the covariates
Computes recursively the denominator of the fixed-effect logit
likelihood for one unit with T
observations and t
ones, as well
as its first derivative wrt the parameter vector beta. The matrix
U
is a vector containing the values of the index function for each
time period and X
is the corresponding matrix of covariates, so
that U[i] = exp(X[i,]*beta).
The first element of the output vector is the likelihood denominator (B in Gail et al's notation) and the subsequent k elements contain the derivative of B wrt beta.
The algorithm is described in Gail et al. (1981), "Likelihood Calculations for Matched Case-Control Studies and Survival Studies with Tied Death Times", Biometrika, Vol. 68 (3), pp. 703-707
DATASET * matrix_dset_plus_lists (const gretl_matrix *m1
,const gretl_matrix *m2
,int **plist1
,int **plist2
,int *err
);