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_var_del (g02df)

Purpose

nag_correg_linregm_var_del (g02df) deletes an independent variable from a general linear regression model.

Syntax

[q, rss, ifail] = g02df(q, indx, rss, 'ip', ip)
[q, rss, ifail] = nag_correg_linregm_var_del(q, indx, rss, 'ip', ip)
Note: the interface to this routine has changed since earlier releases of the toolbox:
Mark 22: ip has been made optional
.

Description

When selecting a linear regression model it is sometimes useful to drop independent variables from the model and to examine the resulting sub-model. nag_correg_linregm_var_del (g02df) updates the QRQR decomposition used in the computation of the linear regression model. The QRQR decomposition may come from nag_correg_linregm_fit (g02da) or nag_correg_linregm_var_add (g02de), or a previous call to nag_correg_linregm_var_del (g02df).
For the general linear regression model with pp independent variables fitted nag_correg_linregm_fit (g02da) or nag_correg_linregm_var_add (g02de) compute a QRQR decomposition of the (weighted) independent variables and form an upper triangular matrix RR and a vector cc. To remove an independent variable RR and cc have to be updated. The column of RR corresponding to the variable to be dropped is removed and the matrix is then restored to upper triangular form by applying a series of Givens rotations. The rotations are then applied to cc. Note only the first pp elements of cc are affected.
The method used means that while the updated values of RR and cc are computed an updated value of QQ from the QRQR decomposition is not available so a call to nag_correg_linregm_var_add (g02de) cannot be made after a call to nag_correg_linregm_var_del (g02df).
nag_correg_linregm_update (g02dd) can be used to calculate the parameter estimates, β̂β^, from the information provided by nag_correg_linregm_var_del (g02df).

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

Parameters

Compulsory Input Parameters

1:     q(ldq,ip + 1ip+1) – double array
ldq, the first dimension of the array, must satisfy the constraint ldqipldqip.
The results of the QRQR decomposition as returned by functions nag_correg_linregm_fit (g02da), nag_correg_linregm_obs_edit (g02dc), nag_correg_linregm_var_add (g02de) or nag_correg_linregm_fit_onestep (g02ee), or previous calls to nag_correg_linregm_var_del (g02df).
2:     indx – int64int32nag_int scalar
Indicates which independent variable is to be deleted from the model.
Constraint: 1indxip1indxip.
3:     rss – double scalar
The residual sum of squares for the full regression.
Constraint: rss0.0rss0.0.

Optional Input Parameters

1:     ip – int64int32nag_int scalar
Default: The dimension of the array q and the first dimension of the array q. (An error is raised if these dimensions are not equal.)
pp, the number of independent variables already in the model.
Constraint: ip1ip1.

Input Parameters Omitted from the MATLAB Interface

ldq wk

Output Parameters

1:     q(ldq,ip + 1ip+1) – double array
ldqipldqip.
The updated QRQR decomposition.
2:     rss – double scalar
The residual sum of squares with the (indx)th variable removed. Note that the residual sum of squares will only be valid if the regression is of full rank, otherwise the residual sum of squares should be obtained using nag_correg_linregm_update (g02dd).
3:     ifail – int64int32nag_int scalar
ifail = 0ifail=0 unless the function detects an error (see [Error Indicators and Warnings]).

Error Indicators and Warnings

Errors or warnings detected by the function:
  ifail = 1ifail=1
On entry,ip < 1ip<1,
orldq < ipldq<ip,
orindx < 1indx<1,
orindx > ipindx>ip,
orrss < 0.0rss<0.0.
  ifail = 2ifail=2
On entry,a diagonal element of RR is zero.

Accuracy

There will inevitably be some loss in accuracy in fitting a model by dropping terms from a more complex model rather than fitting it afresh using nag_correg_linregm_fit (g02da).

Further Comments

None.

Example

function nag_correg_linregm_var_del_example
q = [-27.61466337533981, -3.464101615137754, -12.99038105676658, -23.72909606369362, ...
     -1.732050807568877, -9.093266739736608;
     5.678919489958151, 0.2542949013547401, 5.979130371550697, 8.354057679970834, ...
     1.505235618012764, 8.634198753323211;
     -1.9627798449884, 0.2542949013547401, 0.1680751103201103, -3.676917225029723, ...
     0.1291367593463666, 1.878093517739872;
     -0.06450766139237359, 0.2542949013547401, 0.09396625863990909, ...
     0.3650656216718273, 0.8471065054945881, 2.964167523893722;
     -0.2751502181663455, 0.2542949013547401, 0.01985740695970786, ...
     0.4219147717487357, ...
     0.2859917224052277, 0.4467521874922538];
indx = int64(2);
rss = 0.08406649241210093;
[qOut, rssOut, ifail] = ...
    nag_correg_linregm_var_del(q, indx, rss)
 

qOut =

  -27.6147   -3.4641  -23.7291   -1.7321   -9.0933  -12.9904
    5.9884    0.2543    9.1274    1.3257    7.1460    5.4725
    0.2703    0.2543    0.1681    1.1147    5.6308    1.5656
   -0.3461    0.2543    0.0940    0.3651   -2.0715   -1.7873
   -0.3582    0.2543    0.0199    0.4219    0.2860   -0.3948


rssOut =

    0.2124


ifail =

                    0


function g02df_example
q = [-27.61466337533981, -3.464101615137754, -12.99038105676658, ...
     -23.72909606369362, ...
     -1.732050807568877, -9.093266739736608;
     5.678919489958151, 0.2542949013547401, 5.979130371550697, ...
     8.354057679970834, 1.505235618012764, 8.634198753323211;
     -1.9627798449884, 0.2542949013547401, 0.1680751103201103, ...
     -3.676917225029723, 0.1291367593463666, 1.878093517739872;
     -0.06450766139237359, 0.2542949013547401, 0.09396625863990909, ...
      0.3650656216718273, 0.8471065054945881, 2.964167523893722;
     -0.2751502181663455, 0.2542949013547401, 0.01985740695970786, ...
      0.4219147717487357, ...
     0.2859917224052277, 0.4467521874922538];
indx = int64(2);
rss = 0.08406649241210093;
[qOut, rssOut, ifail] = g02df(q, indx, rss)
 

qOut =

  -27.6147   -3.4641  -23.7291   -1.7321   -9.0933  -12.9904
    5.9884    0.2543    9.1274    1.3257    7.1460    5.4725
    0.2703    0.2543    0.1681    1.1147    5.6308    1.5656
   -0.3461    0.2543    0.0940    0.3651   -2.0715   -1.7873
   -0.3582    0.2543    0.0199    0.4219    0.2860   -0.3948


rssOut =

    0.2124


ifail =

                    0



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–2013