PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_correg_linregm_var_add (g02de)
Purpose
nag_correg_linregm_var_add (g02de) adds a new independent variable to a general linear regression model.
Syntax
[
q,
p,
rss,
ifail] = g02de(
ip,
q,
p,
x, 'n',
n, 'wt',
wt, 'tol',
tol)
[
q,
p,
rss,
ifail] = nag_correg_linregm_var_add(
ip,
q,
p,
x, 'n',
n, 'wt',
wt, 'tol',
tol)
Note: the interface to this routine has changed since earlier releases of the toolbox:
At Mark 23: |
weight was removed from the interface; wt was made optional |
Description
A linear regression model may be built up by adding new independent variables to an existing model.
nag_correg_linregm_var_add (g02de) updates the
decomposition used in the computation of the linear regression model. The
decomposition may come from
nag_correg_linregm_fit (g02da) or a previous call to
nag_correg_linregm_var_add (g02de). The general linear regression model is defined by
where |
is a vector of observations on the dependent variable, |
|
is an by matrix of the independent variables of column rank , |
|
is a vector of length of unknown arguments, |
and |
is a vector of length of unknown random errors such that , where is a known diagonal matrix. |
If , the identity matrix, then least squares estimation is used. If , then for a given weight matrix , weighted least squares estimation is used.
The least squares estimates, of the arguments minimize while the weighted least squares estimates, minimize .
The parameter estimates may be found by computing a
decomposition of
(or
in the weighted case), i.e.,
where
and
is a
by
upper triangular matrix and
is an
by
orthogonal matrix.
If
is of full rank, then
is the solution to
where
(or
) and
is the first
elements of
.
If is not of full rank a solution is obtained by means of a singular value decomposition (SVD) of .
To add a new independent variable, , and have to be updated. The matrix is found such that (or ) is upper triangular. The vector is then updated by multiplying by .
The new independent variable is tested to see if it is linearly related to the existing independent variables by checking that at least one of the values , for , is nonzero.
The new parameter estimates,
, can then be obtained by a call to
nag_correg_linregm_update (g02dd).
The function can be used with , in which case and are initialized.
References
Draper N R and Smith H (1985) Applied Regression Analysis (2nd Edition) Wiley
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Hammarling S (1985) The singular value decomposition in multivariate statistics SIGNUM Newsl. 20(3) 2–25
McCullagh P and Nelder J A (1983) Generalized Linear Models Chapman and Hall
Searle S R (1971) Linear Models Wiley
Parameters
Compulsory Input Parameters
- 1:
– int64int32nag_int scalar
-
, the number of independent variables already in the model.
Constraint:
and .
- 2:
– double array
-
ldq, the first dimension of the array, must satisfy the constraint
.
If
,
q must contain the results of the
decomposition for the model with
arguments as returned by
nag_correg_linregm_fit (g02da) or a previous call to
nag_correg_linregm_var_add (g02de).
If
, the first column of
q should contain the
values of the dependent variable,
.
- 3:
– double array
-
Contains further details of the
decomposition used. The first
ip elements of
p must contain the zeta values for the
decomposition (see
nag_lapack_dgeqrf (f08ae) for details).
The first
ip elements of array
p are provided by
nag_correg_linregm_fit (g02da) or by previous calls to
nag_correg_linregm_var_add (g02de).
- 4:
– double array
-
, the new independent variable.
Optional Input Parameters
- 1:
– int64int32nag_int scalar
-
Default:
the dimension of the array
x and the first dimension of the array
q. (An error is raised if these dimensions are not equal.)
, the number of observations.
Constraint:
.
- 2:
– double array
-
The dimension of the array
wt
must be at least
if
, and at least
otherwise
If provided,
wt must contain the weights to be used.
If , the th observation is not included in the model, in which case the effective number of observations is the number of observations with nonzero weights.
If
wt is not provided the effective number of observations is
.
Constraint:
if , , for .
- 3:
– double scalar
Default:
The value of
tol is used to decide if the new independent variable is linearly related to independent variables already included in the model. If the new variable is linearly related then
is not updated. The smaller the value of
tol the stricter the criterion for deciding if there is a linear relationship.
Constraint:
.
Output Parameters
- 1:
– double array
-
The results of the
decomposition for the model with
arguments:
- the first column of q contains the updated value of ;
- the columns to are unchanged;
- the first elements of column contain the new column of , while the remaining elements contain details of the matrix .
- 2:
– double array
-
The first
ip elements of
p are unchanged and the
th element contains the zeta value for
.
-
The residual sum of squares for the new fitted model.
Note: this will only be valid if the model is of full rank, see
Further Comments.
- 4:
– int64int32nag_int scalar
unless the function detects an error (see
Error Indicators and Warnings).
Error Indicators and Warnings
Note: nag_correg_linregm_var_add (g02de) may return useful information for one or more of the following detected errors or warnings.
Errors or warnings detected by the function:
Cases prefixed with W are classified as warnings and
do not generate an error of type NAG:error_n. See nag_issue_warnings.
-
-
On entry, | , |
or | , |
or | , |
or | , |
or | , |
or | or . |
-
-
On entry, |
and a value of . |
- W
-
The new independent variable is a linear combination of existing variables. The
th column of
q will therefore be null.
-
An unexpected error has been triggered by this routine. Please
contact
NAG.
-
Your licence key may have expired or may not have been installed correctly.
-
Dynamic memory allocation failed.
Accuracy
The accuracy is closely related to the accuracy of
nag_lapack_dormqr (f08ag) which should be consulted for further details.
Further Comments
It should be noted that the residual sum of squares produced by
nag_correg_linregm_var_add (g02de) may not be correct if the model to which the new independent variable is added is not of full rank. In such a case
nag_correg_linregm_update (g02dd) should be used to calculate the residual sum of squares.
Example
A dataset consisting of
observations is read in. The four independent variables are stored in the array
x while the dependent variable is read into the first column of
q. If the character variable
indicates that a mean should be included in the model a variable taking the value
for all observations is set up and fitted. Subsequently, one variable at a time is selected to enter the model as indicated by the input value of
. After the variable has been added the parameter estimates are calculated by
nag_correg_linregm_update (g02dd) and the results printed. This is repeated until the input value of
is
.
Open in the MATLAB editor:
g02de_example
function g02de_example
fprintf('g02de example results\n\n');
x = [1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0 1.0;
1.0 1.5 2.0 2.5 3.0 3.5 4.0 4.5 5.0 5.5 6.0 6.5;
0.0 0.0 0.0 0.0 0.0 0.0 1.0 1.0 1.0 1.0 1.0 1.0;
0.0 0.0 0.0 0.0 0.0 0.0 4.0 4.5 5.0 5.5 6.0 6.5;
1.4 2.2 4.5 6.1 7.1 7.7 8.3 8.6 8.8 9.0 9.3 9.2];
[m,n] = size(x);
y = [4.32; 5.21; 6.49; 7.10; 7.94; 8.53;
8.84; 9.02; 9.27; 9.43; 9.68; 9.83];
q = zeros(n,m+1);
q(:,1) = y;
p = zeros(m*(m+2),1);;
ip = int64(0);
for j = 1:m
[q, p, rss, ifail] = g02de( ...
ip, q, p, x(j,1:n));
ip = ip + 1;
fprintf('\nVariable %4d added\n',ip);
rsst = 0;
[rsst, idf, b, se, covar, svd, irank, p2, ifail] = ...
g02dd(int64(n), ip, q, rsst);
if svd
fprintf('Model not of full rank\n\n');
end
fprintf('Residual sum of squares = %12.4e\n', rsst);
fprintf('Degrees of freedom = %4d\n', idf);
fprintf('\nVariable Parameter estimate Standard error\n\n');
ivar = double([1:ip]');
fprintf('%6d%20.4e%20.4e\n',[ivar b se]');
end
g02de example results
Variable 1 added
Residual sum of squares = 3.6267e+01
Degrees of freedom = 11
Variable Parameter estimate Standard error
1 7.9717e+00 5.2416e-01
Variable 2 added
Residual sum of squares = 4.0164e+00
Degrees of freedom = 10
Variable Parameter estimate Standard error
1 4.4100e+00 4.3756e-01
2 9.4979e-01 1.0599e-01
Variable 3 added
Residual sum of squares = 3.8872e+00
Degrees of freedom = 9
Variable Parameter estimate Standard error
1 4.2236e+00 5.6734e-01
2 1.0554e+00 2.2217e-01
3 -4.1962e-01 7.6695e-01
Variable 4 added
Residual sum of squares = 1.8702e-01
Degrees of freedom = 8
Variable Parameter estimate Standard error
1 2.7605e+00 1.7592e-01
2 1.7057e+00 7.3100e-02
3 4.4575e+00 4.2676e-01
4 -1.3006e+00 1.0338e-01
Variable 5 added
Residual sum of squares = 8.4066e-02
Degrees of freedom = 7
Variable Parameter estimate Standard error
1 3.1440e+00 1.8181e-01
2 9.0748e-01 2.7761e-01
3 2.0790e+00 8.6804e-01
4 -6.1589e-01 2.4530e-01
5 2.9224e-01 9.9810e-02
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015