hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_correg_linregm_update (g02dd)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_correg_linregm_update (g02dd) calculates the regression arguments for a general linear regression model. It is intended to be called after nag_correg_linregm_obs_edit (g02dc), nag_correg_linregm_var_add (g02de) or nag_correg_linregm_var_del (g02df).

Syntax

[rss, idf, b, se, covar, svd, irank, p, ifail] = g02dd(n, ip, q, rss, 'tol', tol)
[rss, idf, b, se, covar, svd, irank, p, ifail] = nag_correg_linregm_update(n, ip, q, rss, 'tol', tol)
Note: the interface to this routine has changed since earlier releases of the toolbox:
At Mark 23: n was made a compulsory input parameter; p was made output only
At Mark 22: n was made optional

Description

A general linear regression model fitted by nag_correg_linregm_fit (g02da) may be adjusted by adding or deleting an observation using nag_correg_linregm_obs_edit (g02dc), adding a new independent variable using nag_correg_linregm_var_add (g02de) or deleting an existing independent variable using nag_correg_linregm_var_del (g02df). Alternatively a model may be constructed by a forward selection procedure using nag_correg_linregm_fit_onestep (g02ee). These functions compute the vector c and the upper triangular matrix R. nag_correg_linregm_update (g02dd) takes these basic results and computes the regression coefficients, β^, their standard errors and their variance-covariance matrix.
If R is of full rank, then β^ is the solution to
Rβ^=c1,  
where c1 is the first p elements of c.
If R is not of full rank a solution is obtained by means of a singular value decomposition (SVD) of R,
R=Q* D 0 0 0 PT,  
where D is a k by k diagonal matrix with nonzero diagonal elements, k being the rank of R, and Q* and P are p by p orthogonal matrices. This gives the solution
β^=P1D-1Q*1Tc1.  
P1 being the first k columns of P, i.e., P=P1P0, and Q*1 being the first k columns of Q*.
Details of the SVD are made available in the form of the matrix P*:
P*= D-1 P1T P0T .  
This will be only one of the possible solutions. Other estimates may be obtained by applying constraints to the arguments. These solutions can be obtained by calling nag_correg_linregm_constrain (g02dk) after calling nag_correg_linregm_update (g02dd). Only certain linear combinations of the arguments will have unique estimates; these are known as estimable functions. These can be estimated using nag_correg_linregm_estfunc (g02dn).
The residual sum of squares required to calculate the standard errors and the variance-covariance matrix can either be input or can be calculated if additional information on c for the whole sample is provided.

References

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
Searle S R (1971) Linear Models Wiley

Parameters

Compulsory Input Parameters

1:     n int64int32nag_int scalar
The number of observations.
Constraint: n1.
2:     ip int64int32nag_int scalar
p, the number of terms in the regression model.
Constraint: ip1.
3:     qldqip+1 – double array
ldq, the first dimension of the array, must satisfy the constraint
  • if rss0.0, ldqn;
  • otherwise ldqip.
4:     rss – double scalar
Either the residual sum of squares or a value less than or equal to 0.0 to indicate that the residual sum of squares is to be calculated by the function.

Optional Input Parameters

1:     tol – double scalar
Default: 0.000001
The value of tol is used to decide if the independent variables are of full rank and, if not, what is the rank of the independent variables. The smaller the value of tol the stricter the criterion for selecting the singular value decomposition. If tol=0.0, the singular value decomposition will never be used, this may cause run time errors or inaccuracies if the independent variables are not of full rank.
Constraint: tol0.0.

Output Parameters

1:     rss – double scalar
If rss0.0 on entry, then on exit rss will contain the residual sum of squares as calculated by nag_correg_linregm_update (g02dd).
If rss was positive on entry, it will be unchanged.
2:     idf int64int32nag_int scalar
The degrees of freedom associated with the residual sum of squares.
3:     bip – double array
The estimates of the p parameters, β^.
4:     seip – double array
The standard errors of the p parameters given in b.
5:     covarip×ip+1/2 – double array
The upper triangular part of the variance-covariance matrix of the p parameter estimates given in b. They are stored packed by column, i.e., the covariance between the parameter estimate given in bi and the parameter estimate given in bj, ji, is stored in covarj×j-1/2+i.
6:     svd – logical scalar
If a singular value decomposition has been performed, svd=true, otherwise svd=false.
7:     irank int64int32nag_int scalar
The rank of the independent variables.
If svd=false, irank=ip.
If svd=true, irank is an estimate of the rank of the independent variables.
irank is calculated as the number of singular values greater than tol× (largest singular value). It is possible for the SVD to be carried out but irank to be returned as ip.
8:     pip×ip+2×ip – double array
Contains details of the singular value decomposition if used.
If svd=false, p is not referenced.
If svd=true, the first ip elements of p will not be referenced, the next ip values contain the singular values. The following ip×ip values contain the matrix P* stored by columns.
9:     ifail int64int32nag_int scalar
ifail=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

Errors or warnings detected by the function:
   ifail=1
On entry,n<1,
orip<1,
orldq<ip,
orldq<n,
ortol<0.0.
   ifail=2
The degrees of freedom for error are less than or equal to 0. In this case the estimates of β are returned but not the standard errors or covariances.
   ifail=3
The singular value decomposition, if used, has failed to converge, see nag_eigen_real_triang_svd (f02wu). This is an unlikely error exit.
   ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
   ifail=-399
Your licence key may have expired or may not have been installed correctly.
   ifail=-999
Dynamic memory allocation failed.

Accuracy

The accuracy of the results will depend on the accuracy of the input R matrix, which may lose accuracy if a large number of observations or variables have been dropped.

Further Comments

None.

Example

A dataset consisting of 12 observations and four independent variables is input and a regression model fitted by calls to nag_correg_linregm_var_add (g02de). The arguments are then calculated by nag_correg_linregm_update (g02dd) and the results printed.
function g02dd_example


fprintf('g02dd example results\n\n');

x = [1, 0, 0, 0;
     0, 0, 0, 1;
     0, 1, 0, 0;
     0, 0, 1, 0;
     0, 0, 0, 1;
     0, 1, 0, 0;
     0, 0, 0, 1;
     1, 0, 0, 0;
     0, 0, 1, 0;
     1, 0, 0, 0;
     0, 0, 1, 0;
     0, 1, 0, 0];
y = [33.63;     39.62;     38.18;     41.46;     38.02;     35.83;
     35.99;     36.58;     42.92;     37.80;     40.43;     37.89];

[n,m]  = size(x);
isx    = ones(m,1,'int64');
ip     = int64(0);
q      = zeros(n,m+1);
lp     = (m*(m+1))/2;
p      = zeros(lp,1);
q(:,1) = y;

% Fit general linear regression model, adding each variable in turn
for j = 1:m
   [q, p, rss, ifail] = g02de( ...
                               ip, q, p, x(:,j));
   ip = ip + 1;
end

% Calculate RSS
rss = 0;
[rss, idf, b, se, covar, svd, irank, p, ifail] = ...
  g02dd(int64(n), ip, q, rss);

% Display results
if svd
  fprintf('Model not of full rank, rank = %4d\n\n', irank);
end
fprintf('Residual sum of squares = %12.4e\n', rss);
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]');


g02dd example results

Residual sum of squares =   2.2227e+01
Degrees of freedom      =    8

Variable   Parameter estimate   Standard error

     1          3.6003e+01          9.6235e-01
     2          3.7300e+01          9.6235e-01
     3          4.1603e+01          9.6235e-01
     4          3.7877e+01          9.6235e-01

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015