arma

   more...
Computes least-square estimates of parameters for an ARMA model.
Synopsis
#include <imsls.h>
float *imsls_f_arma (int n_observations, float z[], int p, int q, ..., 0)
The type double function is imsls_d_arma.
Required Arguments
int n_observations (Input)
Number of observations.
float z[] (Input)
Array of length n_observations containing the observations.
int p (Input)
Number of autoregressive parameters.
int q (Input)
Number of moving average parameters.
Return Value
Pointer to an array of length 1 + p + q with the estimated constant, AR, and MA parameters. If IMSLS_NO_CONSTANT is specified, the 0-th element of this array is 0.0.
Synopsis with Optional Arguments
#include <imsls.h>
float *imsls_f_arma (int n_observations, float z[], int p, int q,
IMSLS_NO_CONSTANT, or
IMSLS_CONSTANT,
IMSLS_AR_LAGS, int ar_lags[],
IMSLS_MA_LAGS, int ma_lags[],
IMSLS_METHOD_OF_MOMENTS, or
IMSLS_LEAST_SQUARES,
IMSLS_BACKCASTING, int maxbc, float  tolerance,
IMSLS_RELATIVE_ERROR, float relative_error,
IMSLS_ABS_FCN_TOL, float afcntol,
IMSLS_GRAD_TOL, float grad_tol,
IMSLS_STEP_TOL, float step_tol,
IMSLS_MAX_ITERATIONS, int max_iterations,
IMSLS_MEAN_ESTIMATE, float *z_mean,
IMSLS_INITIAL_ESTIMATES, float ar[], float ma[],
IMSLS_RESIDUAL, float **residual,
IMSLS_RESIDUAL_USER, float residual[],
IMSLS_PARAM_EST_COV, float **param_est_cov,
IMSLS_PARAM_EST_COV_USER, float param_est_cov[],
IMSLS_AUTOCOV, float **autocov,
IMSLS_AUTOCOV_USER, float autocov[],
IMSLS_SS_RESIDUAL, float *ss_residual,
IMSLS_VAR_NOISE, float *avar,
IMSLS_RETURN_USER, float *constant, float ar[], float ma[],
IMSLS_ARMA_INFO, Imsls_f_arma **arma_info,
0)
Optional Arguments
IMSLS_NO_CONSTANT
or
IMSLS_CONSTANT
If IMSLS_NO_CONSTANT is specified, the time series is not centered about its mean, z_mean. If IMSLS_CONSTANT, the default, is specified, the time series is centered about its mean.
IMSLS_AR_LAGS, int ar_lags[] (Input)
Array of length p containing the order of the autoregressive parameters. The elements of ar_lags must be greater than or equal to 1.
Default: ar_lags = [1, 2, ..., p]
IMSLS_MA_LAGS, int ma_lags[] (Input)
Array of length q containing the order of the moving average parameters. The ma_lags elements must be greater than or equal to 1.
Default: ma_lags = [1, 2, ..., q]
IMSLS_METHOD_OF_MOMENTS
or
IMSLS_LEAST_SQUARES
If IMSLS_METHOD_OF_MOMENTS is specified, the autoregressive and moving average parameters are estimated by a method of moments procedure. If IMSLS_LEAST_SQUARES is specified, the autoregressive and moving average parameters are estimated by a least-squares procedure.
Default: IMSLS_METHOD_OF_MOMENTS is used.
IMSLS_BACKCASTING, int maxbc, float tolerance (Input)
If IMSLS_BACKCASTING is specified, maxbc is the maximum length of backcasting and must be greater than or equal to 0. Argument tolerance is the tolerance level used to determine convergence of the backcast algorithm. Typically, tolerance is set to a fraction of an estimate of the standard deviation of the time series.
Default: maxbc = 10; tolerance = 0.01 × standard deviation of z.
IMSLS_RELATIVE_ERROR, float relative_error (Input)
Stopping criterion for use in the nonlinear equation solver used in both the method of moments and least-squares algorithms.
Default: relative_error = 100 × imsls_f_machine(4).
See documentation for function imsls_f_machine (Chapter 15, Utilities).
IMSLS_ABS_FCN_TOL, float afcntol (Input)
The absolute function tolerance used by the nonlinear least-squares solver that determines the MA parameters in the method of moments. This variable is only needed for non-standard ARMA models.
Default: afcntol = max(10-10,ɛ2/3), where ɛ = imsls_f_machine(4) is the machine precision.
See documentation for function imsls_f_machine (Chapter 15, Utilities).
IMSLS_GRAD_TOL, float grad_tol (Input)
The scaled gradient tolerance used by the nonlinear least-squares solver that determines the MA parameters in the method of moments. This variable is only needed for non-standard ARMA models.
Default: grad_tol = in single, in double, where ɛ = imsls_f_machine(4) is the machine precision.
See documentation for function imsls_f_machine (Chapter 15, Utilities).
IMSLS_STEP_TOL, float step_tol (Input)
The scaled step tolerance used by the nonlinear least-squares solver that determines the MA parameters in the method of moments. This variable is only needed for non-standard ARMA models.
Default: step_tol = ɛ2/3, where ɛ = imsls_f_machine(4) is the machine precision.
See documentation for function imsls_f_machine (Chapter 15, Utilities).
IMSLS_MAX_ITERATIONS, int max_iterations (Input)
Maximum number of iterations allowed in the nonlinear equation solver used in both the method of moments and least-squares algorithms.
Default: max_iterations = 200.
IMSLS_MEAN_ESTIMATE, float *z_mean (Input or Input/Output)
On input, z_mean is an initial estimate of the mean of the time series z. On return, z_mean contains an update of the mean.
If IMSLS_NO_CONSTANT and IMSLS_LEAST_SQUARES are specified, z_mean is not used in parameter estimation.
IMSLS_INITIAL_ESTIMATES, float ar[], float ma[] (Input)
If specified, ar is an array of length p containing preliminary estimates of the autoregressive parameters, and ma is an array of length q containing preliminary estimates of the moving average parameters; otherwise, these are computed internally. IMSLS_INITIAL_ESTIMATES is only applicable if IMSLS_LEAST_SQUARES is also specified.
IMSLS_RESIDUAL, float **residual (Output)
Address of a pointer to an internally allocated array of length
na = (n_observations - max (ar_lags [i]) + maxbc) containing the residuals (including backcasts) at the final parameter estimate point in the first n_observations - max(ar_lags[i]) + nb, where nb is the number of values backcast.
IMSLS_RESIDUAL_USER, float residual[] (Output)
Storage for array residual is provided by the user. See IMSLS_RESIDUAL.
IMSLS_PARAM_EST_COV, float **param_est_cov (Output)
Address of a pointer to an internally allocated array containing the variance-covariance matrix of the estimated ARMA parameters and (optionally) of the estimated mean of series z. The size of the array is np × np, where np = p + q + 1 if z is centered about z_mean, and np = p + q if z is not centered. The ordering of variables in param_est_cov is z_mean, ar, and ma. Argument np must be 1 or larger.
IMSLS_PARAM_EST_COV_USER, float param_est_cov[] (Output)
Storage for array param_est_cov is provided by the user. See IMSLS_PARAM_EST_COV.
IMSLS_AUTOCOV, float **autocov (Output)
Address of a pointer to an array of length p + q + 2 containing the variance and autocovariances of the time series z. Argument autocov[0] contains the variance of the series z. Argument autocov[k] contains the autocovariance at lag k, where k = 0, 1, ..., p + q + 1.
IMSLS_AUTOCOV_USER, float autocov[] (Output)
Storage for array autocov is provided by the user. See IMSLS_AUTOCOV.
IMSLS_SS_RESIDUAL, float *ss_residual (Output)
If specified, ss_residual contains the sum of squares of the random shock, ss_residual = residual[1]2 + ... + residual[na]2, where na is equal to the number of residuals.
IMSLS_VAR_NOISE, float *avar (Output)
If specified, avar contains the innovation variance of the series.
IMSLS_RETURN_USER, float *constant, float ar[], float ma[] (Output)
If specified, constant is the constant parameter estimate, ar is an array of length p containing the final autoregressive parameter estimates, and ma is an array of length q containing the final moving average parameter estimates. If p or q equals zero, a NULL array may be used.
IMSLS_ARMA_INFO, Imsls_f_arma **arma_info (Output)
Address of a pointer to an internally allocated structure of type Imsls_f_arma that contains information necessary in the call to imsls_f_arma_forecast.
Description
Function imsls_f_arma computes estimates of parameters for a nonseasonal ARMA model given a sample of observations, {Wt}, for = 1, 2, ..., n, where n = n_observations. There are two methods, method of moments and least squares, from which to choose. The default is method of moments.
Two methods of parameter estimation, method of moments and least squares, are provided. The user can choose the method of moments algorithm with the optional argument IMSLS_METHOD_OF_MOMENTS. The least-squares algorithm is used if the user specifies IMSLS_LEAST_SQUARES. If the user wishes to use the least-squares algorithm, the preliminary estimates are the method of moments estimates by default. Otherwise, the user can input initial estimates by specifying optional argument IMSLS_INITIAL_ESTIMATES. The following table lists the appropriate optional arguments for both the method of moments and least-squares algorithm:
Method of Moments Only
Least Squares Only
Both Method of Moments and Least Squares
IMSLS_METHOD_OF_MOMENTS
IMSLS_LEAST_SQUARES
IMSLS_RELATIVE_ERROR
IMSLS_ABS_FCN_TOL
IMSLS_CONSTANT
(or IMSLS_NO_CONSTANT)
IMSLS_MAX_ITERATIONS
IMSLS_GRAD_TOL
IMSLS_BACKCASTING
IMSLS_MEAN_ESTIMATE
IMSLS_STEP_TOL
IMSLS_INITIAL_ESTIMATES
IMSLS_AUTOCOV(_USER)
 
IMSLS_RESIDUAL(_USER)
IMSLS_RETURN_USER
 
IMSLS_PARAM_EST_COV(_USER)
IMSLS_ARMA_INFO
 
IMSLS_SS_RESIDUAL
IMSLS_AR_LAGS
 
 
IMSLS_MA_LAGS
 
 
IMSLS_VAR_NOISE
Method of Moments Estimation
The method of moments assumes that the stationary time series {Zt} can be described by a nonseasonal ARMA model of the form
φ(B)Zt = θ0 + θ(B)At, t {0, ±1, ±2, ...}
where B is the backward shift operator, μ is the mean of Zt, and
with p autoregressive and q moving average parameters.
Function imsls_f_arma first orders the AR and MA lags in strictly increasing order. Without loss of generality, it therefore can be assumed that
lφ(1) < lφ(2) < ... < lφ(p), 1 lθ(1) < lθ(2) < ... < lθ(q)
so that the nonseasonal ARMA model is of order (p’, q’), where p’ = lθ (p) and q’ = lθ (q).
In order to keep the notation simple, the following explanations assume the standard ARMA (p, q) model with lφ(i) = i, i = 1,...,p, and lθ (i) = i, i = 1,...,q.
Let μ = z_mean be the estimate of the mean μ of the time series{Zt}, where μ equals the following:
The autocovariance function is estimated by
for k = 0, 1, ..., K, where K = pq. Note that (0) is an estimate of the sample variance.
Given the sample autocovariances, the function computes the method of moments estimates of the autoregressive parameters using the extended Yule-Walker equations as follows:
where
The overall constant θ0 is estimated by the following:
The moving average parameters are estimated based on a system of nonlinear equations given K = p + q + 1 autocovariances, σ(k) for k = 1, ..., K, and p autoregressive parameters φi for i = 1, ..., p.
Let Zʹt = φ(B)Zt. The autocovariances of the derived moving average process Zʹt = θ(B)At are estimated by the following relation:
The iterative procedure for determining the moving average parameters is based on the relation
 
where σ(k) denotes the autocovariance function of the original Zt process.
Let and f =  (f0, f1, ..., fq)T, where
and
The nonlinear system
is solved by a trust-region method. If is the estimate of obtained at the i-th iteration and if a full Newton step is possible, then the new value at the (i + 1)-th iteration is determined by
where
is a square matrix of order q + 1 with entries
The estimation procedure begins with the initial value
and terminates at iteration i when either fi is less than relative_error or i equals max_iterations. The moving average parameter estimates are obtained from the final estimate of by setting
The random shock variance is estimated by the following:
See Box and Jenkins (1976, pp. 498-500) for a description of a function that performs similar computations.
Least-squares Estimation
Suppose the time series {Zt} is generated by a nonseasonal ARMA model of the form,
φ(B) (Zt - μ) = θ(B)At for t  {0, ±1, ±2, …}
where B is the backward shift operator, μ is the mean of Zt, and
with p autoregressive and q moving average parameters. Without loss of generality, the following is assumed:
 lφ (1)  lφ (2)  …  lφ (p)
 lθ (1)  lθ (2)  …  lθ (q)
so that the nonseasonal ARMA model is of order (p’, q’), where p’ = lφ (p) and q’ = lθ (q). Note that the usual hierarchical model assumes the following:
lφ (i) = i, 1  i  p 
lθ (j) = j, 1  j  q 
Consider the sum-of-squares function
where
and T is the backward origin. The random shocks {At} are assumed to be independent and identically distributed
random variables. Hence, the log-likelihood function is given by
where f (μ, φ, θ) is a function of μ, φ, and θ.
For T = 0, the log-likelihood function is conditional on the past values of both Zt and At required to initialize the model. The method of selecting these initial values usually introduces transient bias into the model (Box and Jenkins 1976, pp. 210-211). For T = , this dependency vanishes, and estimation problem concerns maximization of the unconditional log-likelihood function. Box and Jenkins (1976, p. 213) argue that
dominates
The parameter estimates that minimize the sum-of-squares function are called least-squares estimates. For large n, the unconditional least-squares estimates are approximately equal to the maximum likelihood-estimates.
In practice, a finite value of T will enable sufficient approximation of the unconditional sum-of-squares function. The values of [AT] needed to compute the unconditional sum of squares are computed iteratively with initial values of Zt obtained by back forecasting. The residuals (including backcasts), estimate of random shock variance, and covariance matrix of the final parameter estimates also are computed. ARIMA parameters can be computed by using imsls_f_difference with imsls_f_arma.
Examples
Example 1
Consider the Wolfer Sunspot Data (Anderson 1971, p. 660) consisting of the number of sunspots observed each year from 1749 through 1924. The data set for this example consists of the number of sunspots observed from 1770 through 1869. The method of moments estimates
for the ARMA(2, 1) model
where the errors At are independently normally distributed with mean zero and variance .
 
#include <imsls.h>
#include <stdio.h>
 
int main()
{
int p=2, q=1, i, n_observations=100, max_iterations=0;
float w[176][2], z[100], *parameters, relative_error=0.0;
 
imsls_f_data_sets(2, IMSLS_X_COL_DIM,
2, IMSLS_RETURN_USER, w,
0);
for (i=0; i<n_observations; i++) z[i] = w[21+i][1];
 
parameters = imsls_f_arma(n_observations, &z[0], p, q,
IMSLS_RELATIVE_ERROR, relative_error,
IMSLS_MAX_ITERATIONS, max_iterations,
0);
printf("AR estimates are %11.4f and %11.4f.\n",
parameters[1], parameters[2]);
printf("MA estimate is %11.4f.\n", parameters[3]);
}
Output
 
AR estimates are 1.2443 and -0.5751.
MA estimate is -0.1241.
Example 2
The data for this example are the same as that for the initial example. Preliminary method of moments estimates are computed by default, and the method of least squares is used to find the final estimates.
 
#include <imsls.h>
#include <stdio.h>
 
int main()
{
int p=2, q=1, i, n_observations=100;
float w[176][2], z[100], *parameters;
 
imsls_f_data_sets(2, IMSLS_X_COL_DIM,
2, IMSLS_RETURN_USER, w,
0);
for (i=0; i<n_observations; i++) z[i] = w[21+i][1];
 
parameters = imsls_f_arma(n_observations, &z[0], p, q,
IMSLS_LEAST_SQUARES,
0);
printf("AR estimates are %11.4f and %11.4f.\n",
parameters[1], parameters[2]);
printf("MA estimate is %11.4f.\n", parameters[3]);
}
Output
 
AR estimates are 1.5300 and -0.8931.
MA estimate is -0.1324.
Warning Errors
IMSLS_LEAST_SQUARES_FAILED
Least-squares estimation of the parameters has failed to converge. Solution from last iteration is returned. The estimates of the parameters at the last iteration may be used as new starting values.
IMSLS_NEED_POSITIVE_GRADTL
The gradient tolerance must be nonnegative while “grad_tol” = # is given. The algorithm will use “grad_tol” = #.
IMSLS_NEGATIVE_STEP_TOL
The step tolerance must be nonnegative while “step_tol” = # is given. The algorithm will use “step_tol” = #.
IMSLS_NEGATIVE_ABS_FCN_TOL
The absolute function tolerance must be nonnegative while “afcntol” = # is given. The algorithm will use “afcntol” = #.
Fatal Errors
IMSLS_TOO_MANY_CALLS
The number of calls to the function has exceeded “itmax*(“n”+1) = %(i1). The user may try a new initial guess.
IMSLS_INCREASE_ERRREL
The bound for the relative error, “errrel” = %(r1), is too small. No further improvement in the approximate solution is possible. The user should increase “errrel”.
IMSLS_NEW_INITIAL_GUESS
The iteration has not made good progress. The user may try a new initial guess.